Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/semidiscretization/semidiscretization_euler_gravity.jl
5586 views
1
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).
2
# Since these FMAs can increase the performance of many numerical algorithms,
3
# we need to opt-in explicitly.
4
# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.
5
@muladd begin
6
#! format: noindent
7
8
"""
9
ParametersEulerGravity(; background_density=0.0,
10
gravitational_constant=1.0,
11
cfl=1.0,
12
resid_tol=1.0e-4,
13
n_iterations_max=10^4,
14
timestep_gravity=timestep_gravity_erk52_3Sstar!)
15
16
Set up parameters for the gravitational part of a [`SemidiscretizationEulerGravity`](@ref).
17
18
# Arguments
19
- `background_density<:Real`: Constant background/reference density ρ₀ which is subtracted from the (Euler) density
20
in the RHS source term computation of the gravity solver.
21
- `gravitational_constant<:Real`: Gravitational constant G which needs to be in consistent units with the
22
density and velocity fields.
23
- `cfl<:Real`: CFL number used for the pseudo-time stepping to advance the hyperbolic diffusion equations into steady state.
24
- `resid_tol<:Real`: Absolute tolerance for the residual of the hyperbolic diffusion equations which are solved to
25
(approximately) steady state.
26
- `n_iterations_max::Int`: Maximum number of iterations of the pseudo-time gravity solver.
27
If `n_iterations <= 0` the solver will iterate until the residual is less or equal `resid_tol`.
28
This can cause an infinite loop if the solver does not converge!
29
- `timestep_gravity`: Function to advance the gravity solver by one pseudo-time step.
30
There are three optimized methods available:
31
1) `timestep_gravity_erk51_3Sstar!` (first-order),
32
2) `timestep_gravity_erk52_3Sstar!` (second-order),
33
3) `timestep_gravity_erk53_3Sstar!` (third-order).
34
Additionally, `timestep_gravity_carpenter_kennedy_erk54_2N!` (fourth-order) can be used.
35
"""
36
struct ParametersEulerGravity{RealT <: Real, TimestepGravity}
37
background_density :: RealT # aka rho0
38
gravitational_constant :: RealT # aka G
39
cfl :: RealT # CFL number for the gravity solver
40
resid_tol :: RealT # Hyp.-Diff. Eq. steady state tolerance
41
n_iterations_max :: Int # Max. number of iterations of the pseudo-time gravity solver
42
timestep_gravity :: TimestepGravity
43
end
44
45
function ParametersEulerGravity(; background_density = 0.0,
46
gravitational_constant = 1.0,
47
cfl = 1.0,
48
resid_tol = 1.0e-4,
49
n_iterations_max = 10^4,
50
timestep_gravity = timestep_gravity_erk52_3Sstar!)
51
background_density, gravitational_constant, cfl, resid_tol = promote(background_density,
52
gravitational_constant,
53
cfl, resid_tol)
54
return ParametersEulerGravity(background_density, gravitational_constant, cfl,
55
resid_tol,
56
n_iterations_max, timestep_gravity)
57
end
58
59
function Base.show(io::IO, parameters::ParametersEulerGravity)
60
@nospecialize parameters # reduce precompilation time
61
62
print(io, "ParametersEulerGravity(")
63
print(io, "background_density=", parameters.background_density)
64
print(io, ", gravitational_constant=", parameters.gravitational_constant)
65
print(io, ", cfl=", parameters.cfl)
66
print(io, ", n_iterations_max=", parameters.n_iterations_max)
67
print(io, ", timestep_gravity=", parameters.timestep_gravity)
68
print(io, ")")
69
return nothing
70
end
71
function Base.show(io::IO, ::MIME"text/plain", parameters::ParametersEulerGravity)
72
@nospecialize parameters # reduce precompilation time
73
74
if get(io, :compact, false)
75
show(io, parameters)
76
else
77
setup = [
78
"background density (ρ₀)" => parameters.background_density,
79
"gravitational constant (G)" => parameters.gravitational_constant,
80
"CFL (gravity)" => parameters.cfl,
81
"max. #iterations" => parameters.n_iterations_max,
82
"time integrator" => parameters.timestep_gravity
83
]
84
summary_box(io, "ParametersEulerGravity", setup)
85
end
86
end
87
88
"""
89
SemidiscretizationEulerGravity
90
91
A struct containing everything needed to describe a spatial semidiscretization
92
of a the compressible Euler equations with self-gravity, reformulating the
93
Poisson equation for the gravitational potential as steady-state problem of
94
the hyperbolic diffusion equations.
95
- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)
96
"A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics"
97
[arXiv: 2008.10593](https://arXiv.org/abs/2008.10593)
98
"""
99
struct SemidiscretizationEulerGravity{SemiEuler, SemiGravity,
100
Parameters <: ParametersEulerGravity, Cache} <:
101
AbstractSemidiscretization
102
semi_euler :: SemiEuler
103
semi_gravity :: SemiGravity
104
parameters :: Parameters
105
cache :: Cache
106
performance_counter :: PerformanceCounter
107
gravity_counter :: PerformanceCounter
108
end
109
# We assume some properties of the fields of the semidiscretization, e.g.,
110
# the `equations` and the `mesh` should have the same dimension. We check these
111
# properties in the outer constructor defined below. While we could ensure
112
# them even better in an inner constructor, we do not use this approach to
113
# simplify the integration with Adapt.jl for GPU usage, see
114
# https://github.com/trixi-framework/Trixi.jl/pull/2677#issuecomment-3591789921
115
116
"""
117
SemidiscretizationEulerGravity(semi_euler::SemiEuler, semi_gravity::SemiGravity, parameters)
118
119
Construct a semidiscretization of the compressible Euler equations with self-gravity.
120
`parameters` should be given as [`ParametersEulerGravity`](@ref).
121
"""
122
function SemidiscretizationEulerGravity(semi_euler::SemiEuler,
123
semi_gravity::SemiGravity,
124
parameters) where
125
{Mesh,
126
SemiEuler <:
127
SemidiscretizationHyperbolic{Mesh, <:AbstractCompressibleEulerEquations},
128
SemiGravity <:
129
SemidiscretizationHyperbolic{Mesh, <:AbstractHyperbolicDiffusionEquations}}
130
@assert ndims(semi_euler) == ndims(semi_gravity)
131
@assert typeof(semi_euler.mesh) == typeof(semi_gravity.mesh)
132
@assert polydeg(semi_euler.solver) == polydeg(semi_gravity.solver)
133
134
u_ode = compute_coefficients(zero(real(semi_gravity)), semi_gravity)
135
du_ode = similar(u_ode)
136
# Registers for gravity solver, tailored to the 2N and 3S* methods implemented below
137
u_tmp1_ode = similar(u_ode)
138
u_tmp2_ode = similar(u_ode)
139
cache = (; u_ode, du_ode, u_tmp1_ode, u_tmp2_ode)
140
141
performance_counter = PerformanceCounter()
142
gravity_counter = PerformanceCounter()
143
144
return SemidiscretizationEulerGravity{typeof(semi_euler), typeof(semi_gravity),
145
typeof(parameters), typeof(cache)}(semi_euler,
146
semi_gravity,
147
parameters,
148
cache,
149
performance_counter,
150
gravity_counter)
151
end
152
153
function remake(semi::SemidiscretizationEulerGravity;
154
uEltype = real(semi.semi_gravity.solver),
155
semi_euler = semi.semi_euler,
156
semi_gravity = semi.semi_gravity,
157
parameters = semi.parameters)
158
semi_euler = remake(semi_euler, uEltype = uEltype)
159
semi_gravity = remake(semi_gravity, uEltype = uEltype)
160
161
# Recreate cache, i.e., registers for u with e.g. AD datatype
162
u_ode = compute_coefficients(zero(real(semi_gravity)), semi_gravity)
163
du_ode = similar(u_ode)
164
u_tmp1_ode = similar(u_ode)
165
u_tmp2_ode = similar(u_ode)
166
cache = (; u_ode, du_ode, u_tmp1_ode, u_tmp2_ode)
167
168
return SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)
169
end
170
171
function Base.show(io::IO, semi::SemidiscretizationEulerGravity)
172
@nospecialize semi # reduce precompilation time
173
174
print(io, "SemidiscretizationEulerGravity using")
175
print(io, semi.semi_euler)
176
print(io, ", ", semi.semi_gravity)
177
print(io, ", ", semi.parameters)
178
print(io, ", cache(")
179
for (idx, key) in enumerate(keys(semi.cache))
180
idx > 1 && print(io, " ")
181
print(io, key)
182
end
183
print(io, "))")
184
return nothing
185
end
186
187
function Base.show(io::IO, mime::MIME"text/plain", semi::SemidiscretizationEulerGravity)
188
@nospecialize semi # reduce precompilation time
189
190
if get(io, :compact, false)
191
show(io, semi)
192
else
193
summary_header(io, "SemidiscretizationEulerGravity")
194
summary_line(io, "semidiscretization Euler",
195
semi.semi_euler |> typeof |> nameof)
196
show(increment_indent(io), mime, semi.semi_euler)
197
summary_line(io, "semidiscretization gravity",
198
semi.semi_gravity |> typeof |> nameof)
199
show(increment_indent(io), mime, semi.semi_gravity)
200
summary_line(io, "parameters", semi.parameters |> typeof |> nameof)
201
show(increment_indent(io), mime, semi.parameters)
202
summary_footer(io)
203
end
204
end
205
206
# The compressible Euler semidiscretization is considered to be the main semidiscretization.
207
# The hyperbolic diffusion equations part is only used internally to update the gravitational
208
# potential during an rhs! evaluation of the flow solver.
209
@inline function mesh_equations_solver_cache(semi::SemidiscretizationEulerGravity)
210
return mesh_equations_solver_cache(semi.semi_euler)
211
end
212
213
@inline Base.ndims(semi::SemidiscretizationEulerGravity) = ndims(semi.semi_euler)
214
215
@inline Base.real(semi::SemidiscretizationEulerGravity) = real(semi.semi_euler)
216
217
# computes the coefficients of the initial condition
218
@inline function compute_coefficients(t, semi::SemidiscretizationEulerGravity)
219
compute_coefficients!(semi.cache.u_ode, t, semi.semi_gravity)
220
return compute_coefficients(t, semi.semi_euler)
221
end
222
223
# computes the coefficients of the initial condition and stores the Euler part in `u_ode`
224
@inline function compute_coefficients!(u_ode, t, semi::SemidiscretizationEulerGravity)
225
compute_coefficients!(semi.cache.u_ode, t, semi.semi_gravity)
226
return compute_coefficients!(u_ode, t, semi.semi_euler)
227
end
228
229
@inline function calc_error_norms(func, u, t, analyzer,
230
semi::SemidiscretizationEulerGravity, cache_analysis)
231
return calc_error_norms(func, u, t, analyzer, semi.semi_euler, cache_analysis)
232
end
233
234
# Coupled Euler and gravity solver at each Runge-Kutta stage,
235
# corresponding to Algorithm 2 in Schlottke-Lakemper et al. (2020),
236
# https://dx.doi.org/10.1016/j.jcp.2021.110467
237
function rhs!(du_ode, u_ode, semi::SemidiscretizationEulerGravity, t)
238
@unpack semi_euler, semi_gravity, cache = semi
239
240
u_euler = wrap_array(u_ode, semi_euler)
241
du_euler = wrap_array(du_ode, semi_euler)
242
u_gravity = wrap_array(cache.u_ode, semi_gravity)
243
n_elements = size(u_euler)[end]
244
245
time_start = time_ns()
246
247
# standard semidiscretization of the compressible Euler equations
248
@trixi_timeit timer() "Euler solver" rhs!(du_ode, u_ode, semi_euler, t)
249
250
# compute gravitational potential and forces
251
@trixi_timeit timer() "gravity solver" update_gravity!(semi, u_ode)
252
253
# add gravitational source source_terms to the Euler part
254
if ndims(semi_euler) == 1
255
@threaded for element in 1:n_elements
256
@views @. du_euler[2, .., element] -= u_euler[1, .., element] *
257
u_gravity[2, .., element]
258
@views @. du_euler[3, .., element] -= u_euler[2, .., element] *
259
u_gravity[2, .., element]
260
end
261
elseif ndims(semi_euler) == 2
262
@threaded for element in 1:n_elements
263
@views @. du_euler[2, .., element] -= u_euler[1, .., element] *
264
u_gravity[2, .., element]
265
@views @. du_euler[3, .., element] -= u_euler[1, .., element] *
266
u_gravity[3, .., element]
267
@views @. du_euler[4, .., element] -= (u_euler[2, .., element] *
268
u_gravity[2, .., element] +
269
u_euler[3, .., element] *
270
u_gravity[3, .., element])
271
end
272
elseif ndims(semi_euler) == 3
273
@threaded for element in 1:n_elements
274
@views @. du_euler[2, .., element] -= u_euler[1, .., element] *
275
u_gravity[2, .., element]
276
@views @. du_euler[3, .., element] -= u_euler[1, .., element] *
277
u_gravity[3, .., element]
278
@views @. du_euler[4, .., element] -= u_euler[1, .., element] *
279
u_gravity[4, .., element]
280
@views @. du_euler[5, .., element] -= (u_euler[2, .., element] *
281
u_gravity[2, .., element] +
282
u_euler[3, .., element] *
283
u_gravity[3, .., element] +
284
u_euler[4, .., element] *
285
u_gravity[4, .., element])
286
end
287
else
288
error("Number of dimensions $(ndims(semi_euler)) not supported.")
289
end
290
291
runtime = time_ns() - time_start
292
put!(semi.performance_counter, runtime)
293
294
return nothing
295
end
296
297
# TODO: Taal refactor, add some callbacks or so within the gravity update to allow investigating/optimizing it
298
function update_gravity!(semi::SemidiscretizationEulerGravity, u_ode)
299
@unpack semi_euler, semi_gravity, parameters, gravity_counter, cache = semi
300
301
u_euler = wrap_array(u_ode, semi_euler)
302
u_gravity = wrap_array(cache.u_ode, semi_gravity)
303
du_gravity = wrap_array(cache.du_ode, semi_gravity)
304
305
# set up main loop
306
finalstep = false
307
@unpack n_iterations_max, cfl, resid_tol, timestep_gravity = parameters
308
iter = 0
309
tau = zero(real(semi_gravity.solver)) # Pseudo-time
310
311
# iterate gravity solver until convergence or maximum number of iterations are reached
312
@unpack equations = semi_gravity
313
while !finalstep
314
dtau = @trixi_timeit timer() "calculate dtau" begin
315
cfl * max_dt(u_gravity, tau, semi_gravity.mesh,
316
have_constant_speed(equations), equations,
317
semi_gravity.solver, semi_gravity.cache)
318
end
319
320
# evolve solution by one pseudo-time step
321
time_start = time_ns()
322
timestep_gravity(cache, u_euler, tau, dtau, parameters, semi_gravity)
323
runtime = time_ns() - time_start
324
put!(gravity_counter, runtime)
325
326
# update iteration counter
327
iter += 1
328
tau += dtau
329
330
# check if we reached the maximum number of iterations
331
if n_iterations_max > 0 && iter >= n_iterations_max
332
@warn "Max iterations reached: Gravity solver failed to converge!" residual=maximum(abs,
333
@views du_gravity[1,
334
..,
335
:]) tau=tau dtau=dtau
336
finalstep = true
337
end
338
339
# this is an absolute tolerance check
340
if maximum(abs, @views du_gravity[1, .., :]) <= resid_tol
341
finalstep = true
342
end
343
end
344
345
return nothing
346
end
347
348
# Integrate gravity solver for 2N-type low-storage schemes
349
function timestep_gravity_2N!(cache, u_euler, tau, dtau, gravity_parameters,
350
semi_gravity,
351
a, b, c)
352
G = gravity_parameters.gravitational_constant
353
rho0 = gravity_parameters.background_density
354
grav_scale = -4.0 * pi * G
355
356
# Note that `u_ode` is `u_gravity` in `rhs!` above
357
@unpack u_ode, du_ode, u_tmp1_ode = cache
358
n_elements = size(u_euler)[end]
359
360
u_tmp1_ode .= zero(eltype(u_tmp1_ode))
361
du_gravity = wrap_array(du_ode, semi_gravity)
362
363
for stage in eachindex(c)
364
tau_stage = tau + dtau * c[stage]
365
366
# rhs! has the source term for the harmonic problem
367
# We don't need a `@trixi_timeit timer() "rhs!"` here since that's already
368
# included in the `rhs!` call.
369
rhs!(du_ode, u_ode, semi_gravity, tau_stage)
370
371
# Source term: Jeans instability OR coupling convergence test OR blast wave
372
# put in gravity source term proportional to Euler density
373
# OBS! subtract off the background density ρ_0 (spatial mean value)
374
# Note: Adding to `du_gravity` is essentially adding to `du_ode`!
375
@threaded for element in 1:n_elements
376
@views @. du_gravity[1, .., element] += grav_scale *
377
(u_euler[1, .., element] - rho0)
378
end
379
380
a_stage = a[stage]
381
b_stage_dtau = b[stage] * dtau
382
@trixi_timeit timer() "Runge-Kutta step" begin
383
@threaded for idx in eachindex(u_ode)
384
u_tmp1_ode[idx] = du_ode[idx] - u_tmp1_ode[idx] * a_stage
385
u_ode[idx] += u_tmp1_ode[idx] * b_stage_dtau
386
end
387
end
388
end
389
390
return nothing
391
end
392
393
function timestep_gravity_carpenter_kennedy_erk54_2N!(cache, u_euler, tau, dtau,
394
gravity_parameters, semi_gravity)
395
# Coefficients for Carpenter's 5-stage 4th-order low-storage Runge-Kutta method
396
a = SVector(0.0,
397
567301805773.0 / 1357537059087.0,
398
2404267990393.0 / 2016746695238.0,
399
3550918686646.0 / 2091501179385.0,
400
1275806237668.0 / 842570457699.0)
401
b = SVector(1432997174477.0 / 9575080441755.0,
402
5161836677717.0 / 13612068292357.0,
403
1720146321549.0 / 2090206949498.0,
404
3134564353537.0 / 4481467310338.0,
405
2277821191437.0 / 14882151754819.0)
406
c = SVector(0.0,
407
1432997174477.0 / 9575080441755.0,
408
2526269341429.0 / 6820363962896.0,
409
2006345519317.0 / 3224310063776.0,
410
2802321613138.0 / 2924317926251.0)
411
412
return timestep_gravity_2N!(cache, u_euler, tau, dtau, gravity_parameters,
413
semi_gravity, a, b, c)
414
end
415
416
# Integrate gravity solver for 3S*-type low-storage schemes
417
function timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
418
semi_gravity,
419
gamma1, gamma2, gamma3, beta, delta, c)
420
G = gravity_parameters.gravitational_constant
421
rho0 = gravity_parameters.background_density
422
grav_scale = -4 * G * pi
423
424
# Note that `u_ode` is `u_gravity` in `rhs!` above
425
@unpack u_ode, du_ode, u_tmp1_ode, u_tmp2_ode = cache
426
n_elements = size(u_euler)[end]
427
428
u_tmp1_ode .= zero(eltype(u_tmp1_ode))
429
u_tmp2_ode .= u_ode
430
du_gravity = wrap_array(du_ode, semi_gravity)
431
432
for stage in eachindex(c)
433
tau_stage = tau + dtau * c[stage]
434
435
# rhs! has the source term for the harmonic problem
436
# We don't need a `@trixi_timeit timer() "rhs!"` here since that's already
437
# included in the `rhs!` call.
438
rhs!(du_ode, u_ode, semi_gravity, tau_stage)
439
440
# Source term: Jeans instability OR coupling convergence test OR blast wave
441
# put in gravity source term proportional to Euler density
442
# OBS! subtract off the background density ρ_0 around which the Jeans instability is perturbed
443
# Note: Adding to `du_gravity` is essentially adding to `du_ode`!
444
@threaded for element in 1:n_elements
445
@views @. du_gravity[1, .., element] += grav_scale *
446
(u_euler[1, .., element] - rho0)
447
end
448
449
delta_stage = delta[stage]
450
gamma1_stage = gamma1[stage]
451
gamma2_stage = gamma2[stage]
452
gamma3_stage = gamma3[stage]
453
beta_stage_dtau = beta[stage] * dtau
454
@trixi_timeit timer() "Runge-Kutta step" begin
455
@threaded for idx in eachindex(u_ode)
456
# See Algorithm 1 (3S* method) in Schlottke-Lakemper et al. (2020)
457
u_tmp1_ode[idx] += delta_stage * u_ode[idx]
458
u_ode[idx] = (gamma1_stage * u_ode[idx] +
459
gamma2_stage * u_tmp1_ode[idx] +
460
gamma3_stage * u_tmp2_ode[idx] +
461
beta_stage_dtau * du_ode[idx])
462
end
463
end
464
end
465
466
return nothing
467
end
468
469
# First-order, 5-stage, 3S*-storage optimized method
470
function timestep_gravity_erk51_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
471
semi_gravity)
472
# New 3Sstar coefficients optimized for polynomials of degree polydeg=3
473
# and examples/parameters_hypdiff_lax_friedrichs.toml
474
# 5 stages, order 1
475
gamma1 = SVector(0.0000000000000000E+00, 5.2910412316555866E-01,
476
2.8433964362349406E-01, -1.4467571130907027E+00,
477
7.5592215948661057E-02)
478
gamma2 = SVector(1.0000000000000000E+00, 2.6366970460864109E-01,
479
3.7423646095836322E-01, 7.8786901832431289E-01,
480
3.7754129043053775E-01)
481
gamma3 = SVector(0.0000000000000000E+00, 0.0000000000000000E+00,
482
0.0000000000000000E+00, 8.0043329115077388E-01,
483
1.3550099149374278E-01)
484
beta = SVector(1.9189497208340553E-01, 5.4506406707700059E-02,
485
1.2103893164085415E-01, 6.8582252490550921E-01,
486
8.7914657211972225E-01)
487
delta = SVector(1.0000000000000000E+00, 7.8593091509463076E-01,
488
1.2639038717454840E-01, 1.7726945920209813E-01,
489
0.0000000000000000E+00)
490
c = SVector(0.0000000000000000E+00, 1.9189497208340553E-01, 1.9580448818599061E-01,
491
2.4241635859769023E-01, 5.0728347557552977E-01)
492
493
return timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
494
semi_gravity,
495
gamma1, gamma2, gamma3, beta, delta, c)
496
end
497
498
# Second-order, 5-stage, 3S*-storage optimized method
499
function timestep_gravity_erk52_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
500
semi_gravity)
501
# New 3Sstar coefficients optimized for polynomials of degree polydeg=3
502
# and examples/parameters_hypdiff_lax_friedrichs.toml
503
# 5 stages, order 2
504
gamma1 = SVector(0.0000000000000000E+00, 5.2656474556752575E-01,
505
1.0385212774098265E+00, 3.6859755007388034E-01,
506
-6.3350615190506088E-01)
507
gamma2 = SVector(1.0000000000000000E+00, 4.1892580153419307E-01,
508
-2.7595818152587825E-02, 9.1271323651988631E-02,
509
6.8495995159465062E-01)
510
gamma3 = SVector(0.0000000000000000E+00, 0.0000000000000000E+00,
511
0.0000000000000000E+00, 4.1301005663300466E-01,
512
-5.4537881202277507E-03)
513
beta = SVector(4.5158640252832094E-01, 7.5974836561844006E-01,
514
3.7561630338850771E-01, 2.9356700007428856E-02,
515
2.5205285143494666E-01)
516
delta = SVector(1.0000000000000000E+00, 1.3011720142005145E-01,
517
2.6579275844515687E-01, 9.9687218193685878E-01,
518
0.0000000000000000E+00)
519
c = SVector(0.0000000000000000E+00, 4.5158640252832094E-01, 1.0221535725056414E+00,
520
1.4280257701954349E+00, 7.1581334196229851E-01)
521
522
return timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
523
semi_gravity,
524
gamma1, gamma2, gamma3, beta, delta, c)
525
end
526
527
# Third-order, 5-stage, 3S*-storage optimized method
528
function timestep_gravity_erk53_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
529
semi_gravity)
530
# New 3Sstar coefficients optimized for polynomials of degree polydeg=3
531
# and examples/parameters_hypdiff_lax_friedrichs.toml
532
# 5 stages, order 3
533
gamma1 = SVector(0.0000000000000000E+00, 6.9362208054011210E-01,
534
9.1364483229179472E-01, 1.3129305757628569E+00,
535
-1.4615811339132949E+00)
536
gamma2 = SVector(1.0000000000000000E+00, 1.3224582239681788E+00,
537
2.4213162353103135E-01, -3.8532017293685838E-01,
538
1.5603355704723714E+00)
539
gamma3 = SVector(0.0000000000000000E+00, 0.0000000000000000E+00,
540
0.0000000000000000E+00, 3.8306787039991996E-01,
541
-3.5683121201711010E-01)
542
beta = SVector(8.4476964977404881E-02, 3.0834660698015803E-01,
543
3.2131664733089232E-01, 2.8783574345390539E-01,
544
8.2199204703236073E-01)
545
delta = SVector(1.0000000000000000E+00, -7.6832695815481578E-01,
546
1.2497251501714818E-01, 1.4496404749796306E+00,
547
0.0000000000000000E+00)
548
c = SVector(0.0000000000000000E+00, 8.4476964977404881E-02, 2.8110631488732202E-01,
549
5.7093842145029405E-01, 7.2999896418559662E-01)
550
551
return timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,
552
semi_gravity,
553
gamma1, gamma2, gamma3, beta, delta, c)
554
end
555
556
# TODO: Taal decide, where should specific parts like these be?
557
@inline function save_solution_file(u_ode, t, dt, iter,
558
semi::SemidiscretizationEulerGravity,
559
solution_callback,
560
element_variables = Dict{Symbol, Any}();
561
system = "")
562
# If this is called already as part of a multi-system setup (i.e., system is non-empty),
563
# we build a combined system name
564
if !isempty(system)
565
system_euler = system * "_euler"
566
system_gravity = system * "_gravity"
567
else
568
system_euler = "euler"
569
system_gravity = "gravity"
570
end
571
572
u_euler = wrap_array_native(u_ode, semi.semi_euler)
573
filename_euler = save_solution_file(u_euler, t, dt, iter,
574
mesh_equations_solver_cache(semi.semi_euler)...,
575
solution_callback, element_variables,
576
system = system_euler)
577
578
u_gravity = wrap_array_native(semi.cache.u_ode, semi.semi_gravity)
579
filename_gravity = save_solution_file(u_gravity, t, dt, iter,
580
mesh_equations_solver_cache(semi.semi_gravity)...,
581
solution_callback, element_variables,
582
system = system_gravity)
583
584
return filename_euler, filename_gravity
585
end
586
587
@inline function (amr_callback::AMRCallback)(u_ode,
588
semi::SemidiscretizationEulerGravity,
589
t, iter; kwargs...)
590
passive_args = ((semi.cache.u_ode,
591
mesh_equations_solver_cache(semi.semi_gravity)...),)
592
has_changed = amr_callback(u_ode, mesh_equations_solver_cache(semi.semi_euler)...,
593
semi, t, iter;
594
kwargs..., passive_args = passive_args)
595
596
if has_changed
597
new_length = length(semi.cache.u_ode)
598
resize!(semi.cache.du_ode, new_length)
599
resize!(semi.cache.u_tmp1_ode, new_length)
600
resize!(semi.cache.u_tmp2_ode, new_length)
601
end
602
603
return has_changed
604
end
605
end # @muladd
606
607