BVP Functions and Jacobian Types

The SciML ecosystem provides an extensive interface for declaring extra functions associated with the boundary value probems's data. In traditional libraries, there is usually only few options: the Jacobian and the Jacobian of boundary conditions. However, we allow for a large array of pre-computed functions to speed up the calculations. This is offered via the BVPFunction types, which can be passed to the problems.

Function Type Definitions

SciMLBase.BVPFunctionType

DocStringExtensions.TypeDefinition()

A representation of a BVP function f, defined by:

\[\frac{du}{dt} = f(u, p, t)\]

and the constraints:

\[g(u, p, t) = 0\]

If the size of g(u, p, t) is different from the size of u, then the constraints are interpreted as a least squares problem, i.e. the objective function is:

\[\min_{u} \| g_i(u, p, t) \|^2\]

and all of its related functions, such as the Jacobian of f, its gradient with respect to time, and more. For all cases, u0 is the initial condition, p are the parameters, and t is the independent variable.

BVPFunction{iip, specialize}(f, bc;
    mass_matrix = __has_mass_matrix(f) ? f.mass_matrix : I,
    analytic = __has_analytic(f) ? f.analytic : nothing,
    tgrad= __has_tgrad(f) ? f.tgrad : nothing,
    jac = __has_jac(f) ? f.jac : nothing,
    bcjac = __has_jac(bc) ? bc.jac : nothing,
    jvp = __has_jvp(f) ? f.jvp : nothing,
    vjp = __has_vjp(f) ? f.vjp : nothing,
    jac_prototype = __has_jac_prototype(f) ? f.jac_prototype : nothing,
    bcjac_prototype = __has_jac_prototype(bc) ? bc.jac_prototype : nothing,
    sparsity = __has_sparsity(f) ? f.sparsity : jac_prototype,
    paramjac = __has_paramjac(f) ? f.paramjac : nothing,
    syms = nothing,
    indepsym= nothing,
    paramsyms = nothing,
    colorvec = __has_colorvec(f) ? f.colorvec : nothing,
    bccolorvec = __has_colorvec(f) ? bc.colorvec : nothing,
    sys = __has_sys(f) ? f.sys : nothing,
    twopoint::Union{Val, Bool} = Val(false))

Note that both the function f and boundary condition bc are required. f should be given as f(du,u,p,t) or out = f(u,p,t). bc should be given as bc(res, u, p, t). See the section on iip for more details on in-place vs out-of-place handling.

All of the remaining functions are optional for improving or accelerating the usage of f and bc. These include:

  • mass_matrix: the mass matrix M represented in the BVP function. Can be used to determine that the equation is actually a BVP for differential algebraic equation (DAE) if M is singular.
  • analytic(u0,p,t): used to pass an analytical solution function for the analytical solution of the BVP. Generally only used for testing and development of the solvers.
  • tgrad(dT,u,h,p,t) or dT=tgrad(u,p,t): returns $\frac{\partial f(u,p,t)}{\partial t}$
  • jac(J,du,u,p,gamma,t) or J=jac(du,u,p,gamma,t): returns $\frac{df}{du}$
  • bcjac(J,du,u,p,gamma,t) or J=jac(du,u,p,gamma,t): returns $\frac{dbc}{du}$
  • jvp(Jv,v,du,u,p,gamma,t) or Jv=jvp(v,du,u,p,gamma,t): returns the directional derivative $\frac{df}{du} v$
  • vjp(Jv,v,du,u,p,gamma,t) or Jv=vjp(v,du,u,p,gamma,t): returns the adjoint derivative $\frac{df}{du}^\ast v$
  • jac_prototype: a prototype matrix matching the type that matches the Jacobian. For example, if the Jacobian is tridiagonal, then an appropriately sized Tridiagonal matrix can be used as the prototype and integrators will specialize on this structure where possible. Non-structured sparsity patterns should use a SparseMatrixCSC with a correct sparsity pattern for the Jacobian. The default is nothing, which means a dense Jacobian.
  • bcjac_prototype: a prototype matrix matching the type that matches the Jacobian. For example, if the Jacobian is tridiagonal, then an appropriately sized Tridiagonal matrix can be used as the prototype and integrators will specialize on this structure where possible. Non-structured sparsity patterns should use a SparseMatrixCSC with a correct sparsity pattern for the Jacobian. The default is nothing, which means a dense Jacobian.
  • paramjac(pJ,u,p,t): returns the parameter Jacobian $\frac{df}{dp}$.
  • colorvec: a color vector according to the SparseDiffTools.jl definition for the sparsity pattern of the jac_prototype. This specializes the Jacobian construction when using finite differences and automatic differentiation to be computed in an accelerated manner based on the sparsity pattern. Defaults to nothing, which means a color vector will be internally computed on demand when required. The cost of this operation is highly dependent on the sparsity pattern.
  • bccolorvec: a color vector according to the SparseDiffTools.jl definition for the sparsity pattern of the bcjac_prototype. This specializes the Jacobian construction when using finite differences and automatic differentiation to be computed in an accelerated manner based on the sparsity pattern. Defaults to nothing, which means a color vector will be internally computed on demand when required. The cost of this operation is highly dependent on the sparsity pattern.

Additional Options:

  • twopoint: Specify that the BVP is a two-point boundary value problem. Use Val(true) or Val(false) for type stability.

iip: In-Place vs Out-Of-Place

For more details on this argument, see the ODEFunction documentation.

specialize: Controlling Compilation and Specialization

For more details on this argument, see the ODEFunction documentation.

Fields

The fields of the BVPFunction type directly match the names of the inputs.

source
SciMLBase.DynamicalBVPFunctionType

DocStringExtensions.TypeDefinition()

A representation of a dynamical BVP function f, defined by:

\[M \frac{ddu}{dt} = f(du,u,p,t)\]

along with its boundary condition:

\[\]

and all of its related functions, such as the Jacobian of f, its gradient with respect to time, and more. For all cases, u0 is the initial condition, p are the parameters, and t is the independent variable.

Constructor

DynamicalBVPFunction{iip,specialize}(f, bc;
                                    mass_matrix = __has_mass_matrix(f) ? f.mass_matrix : I,
                                    analytic = __has_analytic(f) ? f.analytic : nothing,
                                    tgrad= __has_tgrad(f) ? f.tgrad : nothing,
                                    jac = __has_jac(f) ? f.jac : nothing,
                                    jvp = __has_jvp(f) ? f.jvp : nothing,
                                    vjp = __has_vjp(f) ? f.vjp : nothing,
                                    jac_prototype = __has_jac_prototype(f) ? f.jac_prototype : nothing,
                                    sparsity = __has_sparsity(f) ? f.sparsity : jac_prototype,
                                    paramjac = __has_paramjac(f) ? f.paramjac : nothing,
                                    colorvec = __has_colorvec(f) ? f.colorvec : nothing,
                                    sys = __has_sys(f) ? f.sys : nothing
                                    twopoint::Union{Val, Bool} = Val(false))

Note that only the functions f_i themselves are required. These functions should be given as f_i!(du,du,u,p,t) or ddu = f_i(du,u,p,t). See the section on iip for more details on in-place vs out-of-place handling.

All of the remaining functions are optional for improving or accelerating the usage of f. These include:

  • mass_matrix: the mass matrix M_i represented in the ODE function. Can be used to determine that the equation is actually a differential-algebraic equation (DAE) if M is singular. Note that in this case special solvers are required, see the DAE solver page for more details: https://docs.sciml.ai/DiffEqDocs/stable/solvers/daesolve/. Must be an AbstractArray or an AbstractSciMLOperator. Should be given as a tuple of mass matrices, i.e. `(M1, M_2)` for the mass matrices of equations 1 and 2 respectively.
  • analytic(u0,p,t): used to pass an analytical solution function for the analytical solution of the ODE. Generally only used for testing and development of the solvers.
  • tgrad(dT,du,u,p,t) or dT=tgrad(du,u,p,t): returns $\frac{\partial f(du,u,p,t)}{\partial t}$
  • jac(J,du,u,p,t) or J=jac(du,u,p,t): returns $\frac{df}{du}$
  • jvp(Jv,v,u,p,t) or Jv=jvp(v,u,p,t): returns the directional derivative $\frac{df}{du} v$
  • vjp(Jv,v,u,p,t) or Jv=vjp(v,u,p,t): returns the adjoint derivative $\frac{df}{du}^\ast v$
  • jac_prototype: a prototype matrix matching the type that matches the Jacobian. For example, if the Jacobian is tridiagonal, then an appropriately sized Tridiagonal matrix can be used as the prototype and integrators will specialize on this structure where possible. Non-structured sparsity patterns should use a SparseMatrixCSC with a correct sparsity pattern for the Jacobian. The default is nothing, which means a dense Jacobian.
  • paramjac(pJ,du,u,p,t): returns the parameter Jacobian $\frac{df}{dp}$.
  • colorvec: a color vector according to the SparseDiffTools.jl definition for the sparsity pattern of the jac_prototype. This specializes the Jacobian construction when using finite differences and automatic differentiation to be computed in an accelerated manner based on the sparsity pattern. Defaults to nothing, which means a color vector will be internally computed on demand when required. The cost of this operation is highly dependent on the sparsity pattern.

iip: In-Place vs Out-Of-Place

For more details on this argument, see the ODEFunction documentation.

specialize: Controlling Compilation and Specialization

For more details on this argument, see the ODEFunction documentation.

Fields

The fields of the DynamicalBVPFunction type directly match the names of the inputs.

source