BoundaryValueDiffEqMIRK
Monotonic Implicit Runge Kutta(MIRK) Methods. To only use the MIRK methods form BoundaryVaueDiffEq.jl, you need to install them use the Julia package manager:
using Pkg
Pkg.add("BoundaryValueDiffEqFIRK")
solve(prob::BVProblem, alg, dt; kwargs...)
solve(prob::TwoPointBVProblem, alg, dt; kwargs...)
Full List of Methods
MIRK2
: 2 stage Monotonic Implicit Runge-Kutta method, with defect control adaptivity.MIRK3
: 3 stage Monotonic Implicit Runge-Kutta method, with defect control adaptivity.MIRK4
: 4 stage Monotonic Implicit Runge-Kutta method, with defect control adaptivity.MIRK5
: 5 stage Monotonic Implicit Runge-Kutta method, with defect control adaptivity.MIRK6
: 6 stage Monotonic Implicit Runge-Kutta method, with defect control adaptivity.
Detailed Solvers Explanation
BoundaryValueDiffEqMIRK.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}
}
BoundaryValueDiffEqMIRK.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}
}
BoundaryValueDiffEqMIRK.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}
}
BoundaryValueDiffEqMIRK.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}
}
BoundaryValueDiffEqMIRK.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}
}