-
-
Notifications
You must be signed in to change notification settings - Fork 218
/
Copy pathalgorithms.jl
155 lines (147 loc) · 6.03 KB
/
algorithms.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
@doc generic_solver_docstring(
"""Second order method. Exhibits high stability for real eigenvalues
and is smoothened to allow for moderate sized complex eigenvalues.""",
"ROCK2",
"Stabilized Explicit Method.",
"""Assyr Abdulle, Alexei A. Medovikov. Second Order Chebyshev Methods based on Orthogonal Polynomials.
Numerische Mathematik, 90 (1), pp 1-18, 2001. doi: https://dx.doi.org/10.1007/s002110100292""",
"""
- `min_stages`: The minimum degree of the Chebyshev polynomial.
- `max_stages`: The maximumdegree of the Chebyshev polynomial.
- `eigen_est`: function of the form
`(integrator) -> integrator.eigen_est = upper_bound`,
where `upper_bound` is an estimated upper bound on the spectral radius of the Jacobian matrix.
If `eigen_est` is not provided, `upper_bound` will be estimated using the power iteration.
""",
"""
min_stages = 0,
max_stages = 200,
eigen_est = nothing,
""")
struct ROCK2{E} <: OrdinaryDiffEqAdaptiveAlgorithm
min_stages::Int
max_stages::Int
eigen_est::E
end
function ROCK2(; min_stages = 0, max_stages = 200, eigen_est = nothing)
ROCK2(min_stages, max_stages, eigen_est)
end
@doc generic_solver_docstring(
"""Fourth order method. Exhibits high stability for real eigenvalues
and is smoothened to allow for moderate sized complex eigenvalues.""",
"ROCK4",
"Stabilized Explicit Method.",
"""Assyr Abdulle. Fourth Order Chebyshev Methods With Recurrence Relation. 2002 Society for
Industrial and Applied Mathematics Journal on Scientific Computing, 23(6), pp 2041-2054, 2001.
doi: https://doi.org/10.1137/S1064827500379549""",
"""
- `min_stages`: The minimum degree of the Chebyshev polynomial.
- `max_stages`: The maximumdegree of the Chebyshev polynomial.
- `eigen_est`: function of the form
`(integrator) -> integrator.eigen_est = upper_bound`,
where `upper_bound` is an estimated upper bound on the spectral radius of the Jacobian matrix.
If `eigen_est` is not provided, `upper_bound` will be estimated using the power iteration.
""",
"""
min_stages = 0,
max_stages = 152,
eigen_est = nothing,
""")
struct ROCK4{E} <: OrdinaryDiffEqAdaptiveAlgorithm
min_stages::Int
max_stages::Int
eigen_est::E
end
function ROCK4(; min_stages = 0, max_stages = 152, eigen_est = nothing)
ROCK4(min_stages, max_stages, eigen_est)
end
# SERK methods
for Alg in [:ESERK4, :ESERK5, :RKC]
@eval begin
struct $Alg{E} <: OrdinaryDiffEqAdaptiveAlgorithm
eigen_est::E
end
$Alg(; eigen_est = nothing) = $Alg(eigen_est)
end
end
@doc generic_solver_docstring(
"""Second order method. Exhibits high stability for real eigenvalues.""",
"RKC",
"Stabilized Explicit Method.",
"""B. P. Sommeijer, L. F. Shampine, J. G. Verwer. RKC: An Explicit Solver for Parabolic PDEs,
Journal of Computational and Applied Mathematics, 88(2), pp 315-326, 1998. doi:
https://doi.org/10.1016/S0377-0427(97)00219-7""",
"""
- `eigen_est`: function of the form
`(integrator) -> integrator.eigen_est = upper_bound`,
where `upper_bound` is an estimated upper bound on the spectral radius of the Jacobian matrix.
If `eigen_est` is not provided, `upper_bound` will be estimated using the power iteration.
""",
"""
eigen_est = nothing,
""")
function RKC end
@doc generic_solver_docstring(
"""Fourth order method. Exhibits high stability for real eigenvalues
and is smoothened to allow for moderate sized complex eigenvalues.""",
"ESERK4",
"Stabilized Explicit Method.",
"""J. Martín-Vaquero, B. Kleefeld. Extrapolated stabilized explicit Runge-Kutta methods,
Journal of Computational Physics, 326, pp 141-155, 2016. doi:
https://doi.org/10.1016/j.jcp.2016.08.042.""",
"""
- `eigen_est`: function of the form
`(integrator) -> integrator.eigen_est = upper_bound`,
where `upper_bound` is an estimated upper bound on the spectral radius of the Jacobian matrix.
If `eigen_est` is not provided, `upper_bound` will be estimated using the power iteration.
""",
"""
eigen_est = nothing,
""")
function ESERK4 end
@doc generic_solver_docstring(
"""Fifth order method. Exhibits high stability for real eigenvalues
and is smoothened to allow for moderate sized complex eigenvalues.""",
"ESERK5",
"Stabilized Explicit Method.",
"""J. Martín-Vaquero, A. Kleefeld. ESERK5: A fifth-order extrapolated stabilized explicit Runge-Kutta method,
Journal of Computational and Applied Mathematics, 356, pp 22-36, 2019. doi:
https://doi.org/10.1016/j.cam.2019.01.040.""",
"""
- `eigen_est`: function of the form
`(integrator) -> integrator.eigen_est = upper_bound`,
where `upper_bound` is an estimated upper bound on the spectral radius of the Jacobian matrix.
If `eigen_est` is not provided, `upper_bound` will be estimated using the power iteration.
""",
"""
eigen_est = nothing,
""")
function ESERK5 end
@doc generic_solver_docstring("""Second order method.""",
"SERK2",
"Stabilized Explicit Method.",
"""@article{kleefeld2013serk2v2,
title={SERK2v2: A new second-order stabilized explicit Runge-Kutta method for stiff problems},
author={Kleefeld, B and Martin-Vaquero, J},
journal={Numerical Methods for Partial Differential Equations},
volume={29},
number={1},
pages={170--185},
year={2013},
publisher={Wiley Online Library}}""",
"""
- `controller`: TBD
- `eigen_est`: function of the form
`(integrator) -> integrator.eigen_est = upper_bound`,
where `upper_bound` is an estimated upper bound on the spectral radius of the Jacobian matrix.
If `eigen_est` is not provided, `upper_bound` will be estimated using the power iteration.
""",
"""
controller = :PI
eigen_est = nothing,
""")
struct SERK2{E} <: OrdinaryDiffEqAdaptiveAlgorithm
controller::Symbol
eigen_est::E
end
SERK2(; controller = :PI, eigen_est = nothing) = SERK2(controller, eigen_est)