@@ -190,6 +190,80 @@ function DiffEqBase.ODEFunction{iip}(sys::AbstractODESystem, dvs = states(sys),
190
190
syms = Symbol .(states (sys)))
191
191
end
192
192
193
+ """
194
+ ```julia
195
+ function DiffEqBase.ODEFunctionExpr{iip}(sys::AbstractODESystem, dvs = states(sys),
196
+ ps = parameters(sys);
197
+ version = nothing, tgrad=false,
198
+ jac = false, Wfact = false,
199
+ sparse = false,
200
+ kwargs...) where {iip}
201
+ ```
202
+
203
+ Create a Julia expression for an `ODEFunction` from the [`ODESystem`](@ref).
204
+ The arguments `dvs` and `ps` are used to set the order of the dependent
205
+ variable and parameter vectors, respectively.
206
+ """
207
+ struct ODEFunctionExpr{iip} end
208
+
209
+ function ODEFunctionExpr {iip} (sys:: AbstractODESystem , dvs = states (sys),
210
+ ps = parameters (sys), u0 = nothing ;
211
+ version = nothing , tgrad= false ,
212
+ jac = false , Wfact = false ,
213
+ sparse = false ,linenumbers = false ,
214
+ kwargs... ) where {iip}
215
+
216
+ idx = iip ? 2 : 1
217
+ f = generate_function (sys, dvs, ps; expression= Val{true }, kwargs... )[idx]
218
+ if tgrad
219
+ _tgrad = generate_tgrad (sys, dvs, ps; expression= Val{true }, kwargs... )[idx]
220
+ else
221
+ _tgrad = :nothing
222
+ end
223
+
224
+ if jac
225
+ _jac = generate_jacobian (sys, dvs, ps; sparse = sparse, expression= Val{true }, kwargs... )[idx]
226
+ else
227
+ _jac = :nothing
228
+ end
229
+
230
+ if Wfact
231
+ tmp_Wfact,tmp_Wfact_t = generate_factorized_W (sys, dvs, ps; expression= Val{true }, kwargs... )
232
+ _Wfact = tmp_Wfact[idx]
233
+ _Wfact_t = tmp_Wfact_t[idx]
234
+ else
235
+ _Wfact,_Wfact_t = :nothing ,:nothing
236
+ end
237
+
238
+ M = calculate_massmatrix (sys)
239
+
240
+ _M = (u0 === nothing || M == I) ? M : ArrayInterface. restructure (u0 .* u0' ,M)
241
+
242
+ ex = quote
243
+ f = $ f
244
+ tgrad = $ _tgrad
245
+ jac = $ _jac
246
+ Wfact = $ _Wfact
247
+ Wfact_t = $ _Wfact_t
248
+ M = $ _M
249
+
250
+ ODEFunction {$iip} (f,
251
+ jac = jac,
252
+ tgrad = tgrad,
253
+ Wfact = Wfact,
254
+ Wfact_t = Wfact_t,
255
+ mass_matrix = M,
256
+ syms = $ (Symbol .(states (sys))))
257
+ end
258
+ ! linenumbers ? striplines (ex) : ex
259
+ end
260
+
261
+
262
+ function ODEFunctionExpr (sys:: AbstractODESystem , args... ; kwargs... )
263
+ ODEFunctionExpr {true} (sys, args... ; kwargs... )
264
+ end
265
+
266
+
193
267
function DiffEqBase. ODEProblem (sys:: AbstractODESystem , args... ; kwargs... )
194
268
ODEProblem {true} (sys, args... ; kwargs... )
195
269
end
@@ -225,6 +299,51 @@ function DiffEqBase.ODEProblem{iip}(sys::AbstractODESystem,u0map,tspan,
225
299
ODEProblem {iip} (f,u0,tspan,p;kwargs... )
226
300
end
227
301
302
+ """
303
+ ```julia
304
+ function DiffEqBase.ODEProblemExpr{iip}(sys::AbstractODESystem,u0map,tspan,
305
+ parammap=DiffEqBase.NullParameters();
306
+ version = nothing, tgrad=false,
307
+ jac = false, Wfact = false,
308
+ checkbounds = false, sparse = false,
309
+ linenumbers = true, parallel=SerialForm(),
310
+ kwargs...) where iip
311
+ ```
312
+
313
+ Generates a Julia expression for constructing an ODEProblem from an
314
+ ODESystem and allows for automatically symbolically calculating
315
+ numerical enhancements.
316
+ """
317
+ struct ODEProblemExpr{iip} end
318
+
319
+ function ODEProblemExpr {iip} (sys:: AbstractODESystem ,u0map,tspan,
320
+ parammap= DiffEqBase. NullParameters ();
321
+ version = nothing , tgrad= false ,
322
+ jac = false , Wfact = false ,
323
+ checkbounds = false , sparse = false ,
324
+ linenumbers = false , parallel= SerialForm (),
325
+ kwargs... ) where iip
326
+ dvs = states (sys)
327
+ ps = parameters (sys)
328
+ u0 = varmap_to_vars (u0map,dvs)
329
+ p = varmap_to_vars (parammap,ps)
330
+ f = ODEFunctionExpr {iip} (sys,dvs,ps,u0;tgrad= tgrad,jac= jac,Wfact= Wfact,checkbounds= checkbounds,
331
+ linenumbers= linenumbers,parallel= parallel,
332
+ sparse= sparse)
333
+ ex = quote
334
+ f = $ f
335
+ u0 = $ u0
336
+ tspan = $ tspan
337
+ p = $ p
338
+ ODEProblem (f,u0,tspan,p;$ (kwargs... ))
339
+ end
340
+ ! linenumbers ? striplines (ex) : ex
341
+ end
342
+
343
+ function ODEProblemExpr (sys:: AbstractODESystem , args... ; kwargs... )
344
+ ODEProblemExpr {true} (sys, args... ; kwargs... )
345
+ end
346
+
228
347
229
348
# ## Enables Steady State Problems ###
230
349
function DiffEqBase. SteadyStateProblem (sys:: AbstractODESystem , args... ; kwargs... )
@@ -244,7 +363,7 @@ function DiffEqBase.SteadyStateProblem(sys::AbstractODESystem,u0map,tspan,
244
363
Generates an SteadyStateProblem from an ODESystem and allows for automatically
245
364
symbolically calculating numerical enhancements.
246
365
"""
247
- function DiffEqBase. SteadyStateProblem (sys:: AbstractODESystem ,u0map,
366
+ function DiffEqBase. SteadyStateProblem {iip} (sys:: AbstractODESystem ,u0map,
248
367
parammap= DiffEqBase. NullParameters ();
249
368
version = nothing , tgrad= false ,
250
369
jac = false , Wfact = false ,
@@ -260,3 +379,46 @@ function DiffEqBase.SteadyStateProblem(sys::AbstractODESystem,u0map,
260
379
sparse= sparse)
261
380
SteadyStateProblem (f,u0,p;kwargs... )
262
381
end
382
+
383
+ """
384
+ ```julia
385
+ function DiffEqBase.SteadyStateProblemExpr(sys::AbstractODESystem,u0map,tspan,
386
+ parammap=DiffEqBase.NullParameters();
387
+ version = nothing, tgrad=false,
388
+ jac = false, Wfact = false,
389
+ checkbounds = false, sparse = false,
390
+ linenumbers = true, parallel=SerialForm(),
391
+ kwargs...) where iip
392
+ ```
393
+ Generates a Julia expression for building a SteadyStateProblem from
394
+ an ODESystem and allows for automatically symbolically calculating
395
+ numerical enhancements.
396
+ """
397
+ struct SteadyStateProblemExpr{iip} end
398
+
399
+ function SteadyStateProblemExpr {iip} (sys:: AbstractODESystem ,u0map,
400
+ parammap= DiffEqBase. NullParameters ();
401
+ version = nothing , tgrad= false ,
402
+ jac = false , Wfact = false ,
403
+ checkbounds = false , sparse = false ,
404
+ linenumbers = true , parallel= SerialForm (),
405
+ kwargs... ) where iip
406
+ dvs = states (sys)
407
+ ps = parameters (sys)
408
+ u0 = varmap_to_vars (u0map,dvs)
409
+ p = varmap_to_vars (parammap,ps)
410
+ f = ODEFunctionExpr (sys,dvs,ps,u0;tgrad= tgrad,jac= jac,Wfact= Wfact,checkbounds= checkbounds,
411
+ linenumbers= linenumbers,parallel= parallel,
412
+ sparse= sparse)
413
+ ex = quote
414
+ f = $ f
415
+ u0 = $ u0
416
+ p = $ p
417
+ SteadyStateProblem (f,u0,p;$ (kwargs... ))
418
+ end
419
+ ! linenumbers ? striplines (ex) : ex
420
+ end
421
+
422
+ function SteadyStateProblemExpr (sys:: AbstractODESystem , args... ; kwargs... )
423
+ SteadyStateProblemExpr {true} (sys, args... ; kwargs... )
424
+ end
0 commit comments