Skip to content

Commit 340cf97

Browse files
Anas_AbdelrehimChrisRackauckas
Anas_Abdelrehim
authored andcommitted
fixing cut off args
1 parent d571ca6 commit 340cf97

File tree

1 file changed

+86
-0
lines changed

1 file changed

+86
-0
lines changed

src/common_interface/algorithms.jl

Lines changed: 86 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -23,6 +23,27 @@ The choices for the linear solver are:
2323
:PCG - A preconditioned conjugate gradient method. Only for symmetric linear systems.
2424
:TFQMR - A TFQMR method.
2525
:KLU - A sparse factorization method. Requires that the user specifies a Jacobian. The Jacobian must be set as a sparse matrix in the ODEProblem type.
26+
27+
Example:
28+
29+
CVODE_BDF() # BDF method using Newton + Dense solver
30+
CVODE_BDF(method=:Functional) # BDF method using Functional iterations
31+
CVODE_BDF(linear_solver=:Band,jac_upper=3,jac_lower=3) # Banded solver with nonzero diagonals 3 up and 3 down
32+
CVODE_BDF(linear_solver=:BCG) # Biconjugate gradient method
33+
34+
All of the additional options are available. The full constructor is:
35+
36+
CVODE_BDF(;method=:Newton,linear_solver=:Dense,
37+
jac_upper=0,jac_lower=0,
38+
stored_upper = jac_upper + jac_lower,
39+
non_zero=0,krylov_dim=0,
40+
stability_limit_detect=false,
41+
max_hnil_warns = 10,
42+
max_order = 5,
43+
max_error_test_failures = 7,
44+
max_nonlinear_iters = 3,
45+
max_convergence_failures = 10,
46+
prec = nothing, prec_side = 0)
2647
"""
2748
struct CVODE_BDF{Method, LinearSolver, P, PS} <: SundialsODEAlgorithm{Method, LinearSolver}
2849
jac_upper::Int
@@ -108,6 +129,20 @@ The choices for the linear solver are:
108129
:PCG - A preconditioned conjugate gradient method. Only for symmetric linear systems.
109130
:TFQMR - A TFQMR method.
110131
:KLU - A sparse factorization method. Requires that the user specifies a Jacobian. The Jacobian must be set as a sparse matrix in the ODEProblem type.
132+
133+
All of the additional options are available. The full constructor is:
134+
135+
CVODE_Adams(;method=:Functional,linear_solver=:None,
136+
jac_upper=0,jac_lower=0,
137+
stored_upper = jac_upper + jac_lower,
138+
krylov_dim=0,
139+
stability_limit_detect=false,
140+
max_hnil_warns = 10,
141+
max_order = 12,
142+
max_error_test_failures = 7,
143+
max_nonlinear_iters = 3,
144+
max_convergence_failures = 10,
145+
prec = nothing, psetup = nothing, prec_side = 0)
111146
"""
112147
struct CVODE_Adams{Method, LinearSolver, P, PS} <:
113148
SundialsODEAlgorithm{Method, LinearSolver}
@@ -212,6 +247,34 @@ itable:
212247
ARK436L2SA_DIRK_6_3_4: implicit portion of Kennedy and Carpenter's 4th order method
213248
KVAERNO_7_4_5: Kvaerno's 5th order ESDIRK method
214249
ARK548L2SA_DIRK_8_4_5: implicit portion of Kennedy and Carpenter's 5th order method
250+
251+
These can be set for example via:
252+
253+
ARKODE(Sundials.Explicit(),etable = Sundials.DORMAND_PRINCE_7_4_5)
254+
ARKODE(Sundials.Implicit(),itable = Sundials.KVAERNO_4_2_3)
255+
256+
All of the additional options are available. The full constructor is:
257+
258+
ARKODE(stiffness=Sundials.Implicit();
259+
method=:Newton,linear_solver=:Dense,
260+
jac_upper=0,jac_lower=0,stored_upper = jac_upper+jac_lower,
261+
non_zero=0,krylov_dim=0,
262+
max_hnil_warns = 10,
263+
max_error_test_failures = 7,
264+
max_nonlinear_iters = 3,
265+
max_convergence_failures = 10,
266+
predictor_method = 0,
267+
nonlinear_convergence_coefficient = 0.1,
268+
dense_order = 3,
269+
order = 4,
270+
set_optimal_params = false,
271+
crdown = 0.3,
272+
dgmax = 0.2,
273+
rdiv = 2.3,
274+
msbp = 20,
275+
adaptivity_method = 0,
276+
prec = nothing, psetup = nothing, prec_side = 0
277+
)
215278
"""
216279
struct ARKODE{Method, LinearSolver, MassLinearSolver, T, T1, T2, P, PS} <:
217280
SundialsODEAlgorithm{Method, LinearSolver}
@@ -367,6 +430,29 @@ linearsolver - This is the linear solver which is used in the Newton iterations.
367430
:PCG - A preconditioned conjugate gradient method. Only for symmetric linear systems.
368431
:TFQMR - A TFQMR method.
369432
:KLU - A sparse factorization method. Requires that the user specifies a Jacobian. The Jacobian must be set as a sparse matrix in the ODEProblem type.
433+
434+
Example:
435+
436+
IDA() # Newton + Dense solver
437+
IDA(linear_solver=:Band,jac_upper=3,jac_lower=3) # Banded solver with nonzero diagonals 3 up and 3 down
438+
IDA(linear_solver=:BCG) # Biconjugate gradient method
439+
440+
All of the additional options are available. The constructor is:
441+
442+
IDA(;linear_solver=:Dense,jac_upper=0,jac_lower=0,krylov_dim=0,
443+
max_order = 5,
444+
max_error_test_failures = 7,
445+
max_nonlinear_iters = 3,
446+
nonlinear_convergence_coefficient = 0.33,
447+
nonlinear_convergence_coefficient_ic = 0.0033,
448+
max_num_steps_ic = 5,
449+
max_num_jacs_ic = 4,
450+
max_num_iters_ic = 10,
451+
max_num_backs_ic = 100,
452+
use_linesearch_ic = true,
453+
max_convergence_failures = 10,
454+
init_all = false,
455+
prec = nothing, psetup = nothing, prec_side = 0)
370456
"""
371457
struct IDA{LinearSolver, P, PS} <: SundialsDAEAlgorithm{LinearSolver}
372458
jac_upper::Int

0 commit comments

Comments
 (0)