BoundaryValueDiffEq.BVPM2
— TypeBVPM2(; max_num_subintervals = 3000, method_choice = 4, diagnostic_output = 1,
error_control = 1, singular_term = nothing)
BVPM2(max_num_subintervals::Int, method_choice::Int, diagnostic_output::Int,
error_control::Int, singular_term)
Fortran code for solving two-point boundary value problems. For detailed documentation, see ODEInterface.jl.
Keyword Arguments:
- `max_num_subintervals`: Number of maximal subintervals, default as 3000.
- `method_choice`: Choice for IVP-solvers, default as Runge-Kutta method of order 4,
available choices:
- `2`: Runge-Kutta method of order 2.
- `4`: Runge-Kutta method of order 4.
- `6`: Runge-Kutta method of order 6.
- `diagnostic_output`: Diagnostic output for BVPM2, default as non printout, available
choices:
- `-1`: Full diagnostic printout.
- `0`: Selected printout.
- `1`: No printout.
- `error_control`: Determines the error-estimation for which RTOL is used, default as
defect control, available choices:
- `1`: Defect control.
- `2`: Global error control.
- `3`: Defect and then global error control.
- `4`: Linear combination of defect and global error control.
- `singular_term`: either nothing if the ODEs have no singular terms at the left
boundary or a constant (d,d) matrix for the
singular term.
Only available if the ODEInterface
package is loaded.
BoundaryValueDiffEq.BVPSOL
— TypeBVPSOL(; bvpclass = 2, sol_method = 0, odesolver = nothing)
BVPSOL(bvpclass::Int, sol_methods::Int, odesolver)
A FORTRAN77 code which solves highly nonlinear two point boundary value problems using a local linear solver (condensing algorithm) or a global sparse linear solver for the solution of the arising linear subproblems, by Peter Deuflhard, Georg Bader, Lutz Weimann. For detailed documentation, see ODEInterface.jl.
Keyword Arguments
- `bvpclass`: Boundary value problem classification, default as highly nonlinear with
bad initial data, available choices:
- `0`: Linear boundary value problem.
- `1`: Nonlinear with good initial data.
- `2`: Highly Nonlinear with bad initial data.
- `3`: Highly nonlinear with bad initial data and initial rank reduction to
seperable linear boundary conditions.
- `sol_method`: Switch for solution methods, default as local linear solver with
condensing algorithm, available choices:
- `0`: Use local linear solver with condensing algorithm.
- `1`: Use global sparse linear solver.
- `odesolver`: Either `nothing` or ode-solver(dopri5, dop853, seulex, etc.).
Only available if the ODEInterface
package is loaded.
BoundaryValueDiffEq.COLNEW
— TypeCOLNEW(; bvpclass = 1, collocationpts = 7, diagnostic_output = 1,
max_num_subintervals = 3000)
COLNEW(bvpclass::Int, collocationpts::Int, diagnostic_output::Int,
max_num_subintervals::Int)
Keyword Arguments:
bvpclass
: Boundary value problem classification, default as nonlinear and "extra sensitive", available choices:0
: Linear boundary value problem.1
: Nonlinear and regular.2
: Nonlinear and "extra sensitive" (first relax factor is rstart and the nonlinear iteration does not rely on past convergence).3
: fail-early: return immediately upon: a) two successive non-convergences. b) after obtaining an error estimate for the first time.
collocationpts
: Number of collocation points per subinterval. Require orders[i] ≤ k ≤ 7, default as 7diagnostic_output
: Diagnostic output for COLNEW, default as no printout, available choices:-1
: Full diagnostic printout.0
: Selected printout.1
: No printout.
max_num_subintervals
: Number of maximal subintervals, default as 3000.
A Fortran77 code solves a multi-points boundary value problems for a mixed order system of ODEs. It incorporates a new basis representation replacing b-splines, and improvements for the linear and nonlinear algebraic equation solvers.
Only supports two-point boundary value problems.
Only available if the ODEInterface
package is loaded.
BoundaryValueDiffEq.MIRK2
— TypeMIRK2(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm(),
defect_threshold = 0.1, max_num_subintervals = 3000)
2th order Monotonic Implicit Runge Kutta method.
Keyword Arguments
nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used. Note that any autodiff argument for the solver will be ignored and a custom jacobian algorithm will be used.jac_alg
: Jacobian Algorithm used for the nonlinear solver. Defaults toBVPJacobianAlgorithm()
, which automatically decides the best algorithm to use based on the input types and problem type.- For
TwoPointBVProblem
, onlydiffmode
is used (defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
). - For
BVProblem
,bc_diffmode
andnonbc_diffmode
are used. Fornonbc_diffmode
defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
. Forbc_diffmode
, defaults toAutoForwardDiff
if possible elseAutoFiniteDiff
.
- For
defect_threshold
: Threshold for defect control.max_num_subintervals
: Number of maximal subintervals, default as 3000.
For type-stability, the chunksizes for ForwardDiff ADTypes in BVPJacobianAlgorithm
must be provided.
References
@article{Enright1996RungeKuttaSW,
title={Runge-Kutta Software with Defect Control for Boundary Value ODEs},
author={Wayne H. Enright and Paul H. Muir},
journal={SIAM J. Sci. Comput.},
year={1996},
volume={17},
pages={479-497}
}
BoundaryValueDiffEq.MIRK3
— TypeMIRK3(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm(),
defect_threshold = 0.1, max_num_subintervals = 3000)
3th order Monotonic Implicit Runge Kutta method.
Keyword Arguments
nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used. Note that any autodiff argument for the solver will be ignored and a custom jacobian algorithm will be used.jac_alg
: Jacobian Algorithm used for the nonlinear solver. Defaults toBVPJacobianAlgorithm()
, which automatically decides the best algorithm to use based on the input types and problem type.- For
TwoPointBVProblem
, onlydiffmode
is used (defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
). - For
BVProblem
,bc_diffmode
andnonbc_diffmode
are used. Fornonbc_diffmode
defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
. Forbc_diffmode
, defaults toAutoForwardDiff
if possible elseAutoFiniteDiff
.
- For
defect_threshold
: Threshold for defect control.max_num_subintervals
: Number of maximal subintervals, default as 3000.
For type-stability, the chunksizes for ForwardDiff ADTypes in BVPJacobianAlgorithm
must be provided.
References
@article{Enright1996RungeKuttaSW,
title={Runge-Kutta Software with Defect Control for Boundary Value ODEs},
author={Wayne H. Enright and Paul H. Muir},
journal={SIAM J. Sci. Comput.},
year={1996},
volume={17},
pages={479-497}
}
BoundaryValueDiffEq.MIRK4
— TypeMIRK4(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm(),
defect_threshold = 0.1, max_num_subintervals = 3000)
4th order Monotonic Implicit Runge Kutta method.
Keyword Arguments
nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used. Note that any autodiff argument for the solver will be ignored and a custom jacobian algorithm will be used.jac_alg
: Jacobian Algorithm used for the nonlinear solver. Defaults toBVPJacobianAlgorithm()
, which automatically decides the best algorithm to use based on the input types and problem type.- For
TwoPointBVProblem
, onlydiffmode
is used (defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
). - For
BVProblem
,bc_diffmode
andnonbc_diffmode
are used. Fornonbc_diffmode
defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
. Forbc_diffmode
, defaults toAutoForwardDiff
if possible elseAutoFiniteDiff
.
- For
defect_threshold
: Threshold for defect control.max_num_subintervals
: Number of maximal subintervals, default as 3000.
For type-stability, the chunksizes for ForwardDiff ADTypes in BVPJacobianAlgorithm
must be provided.
References
@article{Enright1996RungeKuttaSW,
title={Runge-Kutta Software with Defect Control for Boundary Value ODEs},
author={Wayne H. Enright and Paul H. Muir},
journal={SIAM J. Sci. Comput.},
year={1996},
volume={17},
pages={479-497}
}
BoundaryValueDiffEq.MIRK5
— TypeMIRK5(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm(),
defect_threshold = 0.1, max_num_subintervals = 3000)
5th order Monotonic Implicit Runge Kutta method.
Keyword Arguments
nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used. Note that any autodiff argument for the solver will be ignored and a custom jacobian algorithm will be used.jac_alg
: Jacobian Algorithm used for the nonlinear solver. Defaults toBVPJacobianAlgorithm()
, which automatically decides the best algorithm to use based on the input types and problem type.- For
TwoPointBVProblem
, onlydiffmode
is used (defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
). - For
BVProblem
,bc_diffmode
andnonbc_diffmode
are used. Fornonbc_diffmode
defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
. Forbc_diffmode
, defaults toAutoForwardDiff
if possible elseAutoFiniteDiff
.
- For
defect_threshold
: Threshold for defect control.max_num_subintervals
: Number of maximal subintervals, default as 3000.
For type-stability, the chunksizes for ForwardDiff ADTypes in BVPJacobianAlgorithm
must be provided.
References
@article{Enright1996RungeKuttaSW,
title={Runge-Kutta Software with Defect Control for Boundary Value ODEs},
author={Wayne H. Enright and Paul H. Muir},
journal={SIAM J. Sci. Comput.},
year={1996},
volume={17},
pages={479-497}
}
BoundaryValueDiffEq.MIRK6
— TypeMIRK6(; nlsolve = NewtonRaphson(), jac_alg = BVPJacobianAlgorithm(),
defect_threshold = 0.1, max_num_subintervals = 3000)
6th order Monotonic Implicit Runge Kutta method.
Keyword Arguments
nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used. Note that any autodiff argument for the solver will be ignored and a custom jacobian algorithm will be used.jac_alg
: Jacobian Algorithm used for the nonlinear solver. Defaults toBVPJacobianAlgorithm()
, which automatically decides the best algorithm to use based on the input types and problem type.- For
TwoPointBVProblem
, onlydiffmode
is used (defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
). - For
BVProblem
,bc_diffmode
andnonbc_diffmode
are used. Fornonbc_diffmode
defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
. Forbc_diffmode
, defaults toAutoForwardDiff
if possible elseAutoFiniteDiff
.
- For
defect_threshold
: Threshold for defect control.max_num_subintervals
: Number of maximal subintervals, default as 3000.
For type-stability, the chunksizes for ForwardDiff ADTypes in BVPJacobianAlgorithm
must be provided.
References
@article{Enright1996RungeKuttaSW,
title={Runge-Kutta Software with Defect Control for Boundary Value ODEs},
author={Wayne H. Enright and Paul H. Muir},
journal={SIAM J. Sci. Comput.},
year={1996},
volume={17},
pages={479-497}
}
BoundaryValueDiffEq.MultipleShooting
— TypeMultipleShooting(; nshoots::Int, ode_alg = nothing, nlsolve = nothing,
grid_coarsening = true, jac_alg = nothing)
MultipleShooting(nshoots::Int; kwargs...)
MultipleShooting(nshoots::Int, ode_alg; kwargs...)
MultipleShooting(nshoots::Int, ode_alg, nlsolve; kwargs...)
Multiple Shooting method, reduces BVP to an initial value problem and solves the IVP. Significantly more stable than Single Shooting.
Arguments
nshoots
: Number of shooting points.ode_alg
: ODE algorithm to use for solving the IVP. Any solver which conforms to the SciMLODEProblem
interface can be used! (Defaults tonothing
which will use poly-algorithm ifDifferentialEquations.jl
is loaded else this must be supplied)nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used.jac_alg
: Jacobian Algorithm used for the nonlinear solver. Defaults toBVPJacobianAlgorithm()
, which automatically decides the best algorithm to use based on the input types and problem type.- For
TwoPointBVProblem
, onlydiffmode
is used (defaults toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
). - For
BVProblem
,bc_diffmode
andnonbc_diffmode
are used. Fornonbc_diffmode
we default toAutoSparse(AutoForwardDiff())
if possible elseAutoSparse(AutoFiniteDiff())
. Forbc_diffmode
, we default toAutoForwardDiff
if possible elseAutoFiniteDiff
.
- For
grid_coarsening
: Coarsening the multiple-shooting grid to generate a stable IVP solution. Possible Choices:true
: Halve the grid size, till we reach a grid size of 1.false
: Do not coarsen the grid. Solve a Multiple Shooting Problem and finally solve a Single Shooting Problem.AbstractVector{<:Int}
orNtuple{N, <:Integer}
: Use the provided grid coarsening. For example, ifnshoots = 10
andgrid_coarsening = [5, 2]
, then the grid will be coarsened to[5, 2]
. Note that1
should not be present in the grid coarsening.Function
: Takes the current number of shooting points and returns the next number of shooting points. For example, ifnshoots = 10
andgrid_coarsening = n -> n ÷ 2
, then the grid will be coarsened to[5, 2]
.
BoundaryValueDiffEq.Shooting
— TypeShooting(ode_alg; kwargs...)
Shooting(ode_alg, nlsolve; kwargs...)
Shooting(; ode_alg = nothing, nlsolve = nothing, jac_alg = nothing)
Single shooting method, reduces BVP to an initial value problem and solves the IVP.
Arguments
ode_alg
: ODE algorithm to use for solving the IVP. Any solver which conforms to the SciMLODEProblem
interface can be used! (Defaults tonothing
which will use poly-algorithm ifDifferentialEquations.jl
is loaded else this must be supplied)nlsolve
: Internal Nonlinear solver. Any solver which conforms to the SciMLNonlinearProblem
interface can be used. Note that any autodiff argument for the solver will be ignored and a custom jacobian algorithm will be used.jac_alg
: Jacobian Algorithm used for the Nonlinear Solver. If this is not set, we check ifnlsolve.ad
exists and is not nothing. If it is, we use that to construct the jacobian. If not, we try to use the best algorithm based on the input types and problem type. IfBVPJacobianAlgorithm
is provided, onlydiffmode
is used (defaults toAutoForwardDiff
if possible elseAutoFiniteDiff
).
BoundaryValueDiffEq.__expand_cache!
— Method__expand_cache!(cache::MIRKCache)
After redistributing or halving the mesh, this function expands the required vectors to match the length of the new mesh.
BoundaryValueDiffEq.__extract_mesh
— Method__extract_mesh(u₀, t₀, t₁, n)
Takes the input initial guess and returns the mesh.
BoundaryValueDiffEq.__extract_u0
— Method__extract_u0(u₀, t₀)
Takes the input initial guess and returns the value at the starting mesh point.
BoundaryValueDiffEq.__flatten_initial_guess
— Method__flatten_initial_guess(u₀) -> Union{AbstractMatrix, AbstractVector, Nothing}
Flattens the initial guess into a matrix. For a function u₀
, it returns nothing
. For no initial guess, it returns vec(u₀)
.
BoundaryValueDiffEq.__generate_sparse_jacobian_prototype
— Method__generate_sparse_jacobian_prototype(::MIRKCache, ya, yb, M, N)
__generate_sparse_jacobian_prototype(::MIRKCache, _, ya, yb, M, N)
__generate_sparse_jacobian_prototype(::MIRKCache, ::TwoPointBVProblem, ya, yb, M, N)
Generate a prototype of the sparse Jacobian matrix for the BVP problem with row and column coloring.
If the problem is a TwoPointBVProblem, then this is the complete Jacobian, else it only computes the sparse part excluding the contributions from the boundary conditions.
BoundaryValueDiffEq.__generate_sparse_jacobian_prototype
— Method__generate_sparse_jacobian_prototype(::MultipleShooting, ::StandardBVProblem,
bcresid_prototype, u0, N::Int, nshoots::Int)
__generate_sparse_jacobian_prototype(::MultipleShooting, ::TwoPointBVProblem,
bcresid_prototype, u0, N::Int, nshoots::Int)
Returns a 3-Tuple:
- Entire Jacobian Prototype (if Two-Point Problem) else
nothing
. - Sparse Non-BC Part Jacobian Prototype along with the column and row color vectors.
- Sparse BC Part Jacobian Prototype along with the column and row color vectors (if Two-Point Problem) else
nothing
.
BoundaryValueDiffEq.__has_initial_guess
— Method__has_initial_guess(u₀) -> Bool
Returns true
if the input has an initial guess.
BoundaryValueDiffEq.__initial_guess_length
— Method__initial_guess_length(u₀) -> Int
Returns the length of the initial guess. If the initial guess is a function or no initial guess is supplied, it returns -1
.
BoundaryValueDiffEq.__initial_guess_on_mesh
— Method__initial_guess_on_mesh(u₀, mesh, p, alias_u0::Bool)
Returns the initial guess on the mesh. For DiffEqArray
assumes that the mesh is the same as the mesh of the DiffEqArray
.
If alias_u0
is set to true
, we try our best to minimize copies. This means that u₀
or parts of it will get mutated.
BoundaryValueDiffEq.concrete_jacobian_algorithm
— Methodconcrete_jacobian_algorithm(jac_alg, prob, alg)
concrete_jacobian_algorithm(jac_alg, problem_type, prob, alg)
If user provided all the required fields, then return the user provided algorithm. Otherwise, based on the problem type and the algorithm, decide the missing fields.
For example, for TwoPointBVProblem
, the bc_diffmode
is set to AutoSparse(AutoForwardDiff())
while for StandardBVProblem
, the bc_diffmode
is set to AutoForwardDiff()
.
BoundaryValueDiffEq.defect_estimate!
— Methoddefect_estimate!(cache::MIRKCache)
defectestimate use the discrete solution approximation Y, plus stages of the RK method in 'kdiscrete', plus some new stages in 'k_interp' to construct an interpolant
BoundaryValueDiffEq.half_mesh!
— Methodhalf_mesh!(mesh, mesh_dt)
half_mesh!(cache::MIRKCache)
The input mesh has length of n + 1
. Divide the original subinterval into two equal length subinterval. The mesh
and mesh_dt
are modified in place.
BoundaryValueDiffEq.interp_eval!
— Methodinterp_eval!(y::AbstractArray, cache::MIRKCache, t)
After we construct an interpolant, we use interp_eval to evaluate it.
BoundaryValueDiffEq.interp_setup!
— Methodinterp_setup!(cache::MIRKCache)
interp_setup!
prepare the extra stages in kiinterp for interpolant construction. Here, the kiinterp is the stages in one subinterval.
BoundaryValueDiffEq.interp_weights
— Functioninterp_weights(τ, alg)
interp_weights: solver-specified interpolation weights and its first derivative
BoundaryValueDiffEq.interval
— Methodinterval(mesh, t)
Find the interval that t
belongs to in mesh
. Assumes that mesh
is sorted.
BoundaryValueDiffEq.mesh_selector!
— Methodmesh_selector!(cache::MIRKCache)
Generate new mesh based on the defect.
BoundaryValueDiffEq.redistribute!
— Methodredistribute!(cache::MIRKCache, Nsub_star, ŝ, mesh, mesh_dt)
Generate a new mesh based on the ŝ
.
BoundaryValueDiffEq.sum_stages!
— Functionsum_stages!(cache::MIRKCache, w, w′, i::Int)
sum_stages add the discrete solution, RK method stages and extra stages to construct interpolant.