Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
trixi-framework
GitHub Repository: trixi-framework/Trixi.jl
Path: blob/main/src/callbacks_step/glm_speed.jl
2055 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
GlmSpeedCallback(; glm_scale=0.5, cfl, semi_indices=())
10
11
Update the divergence cleaning wave speed `c_h` according to the time step
12
computed in [`StepsizeCallback`](@ref) for the ideal GLM-MHD equations, the multi-component
13
GLM-MHD equations, and the multi-ion GLM-MHD equations.
14
The `cfl` number should be set to the same value as for the time step size calculation.
15
As for standard [`StepsizeCallback`](@ref) `cfl` can be either set to a `Real` number or
16
a function of time `t` returning a `Real` number.
17
18
The `glm_scale` ensures that the GLM wave speed is lower than the fastest physical waves in the MHD
19
solution and should thus be set to a value within the interval [0,1]. Note that `glm_scale = 0`
20
deactivates the divergence cleaning.
21
22
In case of coupled semidiscretizations, specify for which `semi_index`, i.e. index of the
23
semidiscretization, the divergence cleaning should be applied. See also
24
[`SemidiscretizationCoupled`](@ref).
25
Note: `SemidiscretizationCoupled` and all related features are considered experimental and
26
may change at any time.
27
"""
28
struct GlmSpeedCallback{RealT <: Real, CflType}
29
glm_scale::RealT
30
cfl::CflType
31
semi_indices::Vector{Int}
32
end
33
34
function Base.show(io::IO, cb::DiscreteCallback{<:Any, <:GlmSpeedCallback})
35
@nospecialize cb # reduce precompilation time
36
37
glm_speed_callback = cb.affect!
38
@unpack glm_scale, cfl, semi_indices = glm_speed_callback
39
print(io, "GlmSpeedCallback(glm_scale=", glm_scale, ", cfl=", cfl, "semi_indices=",
40
semi_indices, ")")
41
end
42
43
function Base.show(io::IO, ::MIME"text/plain",
44
cb::DiscreteCallback{<:Any, <:GlmSpeedCallback})
45
@nospecialize cb # reduce precompilation time
46
47
if get(io, :compact, false)
48
show(io, cb)
49
else
50
glm_speed_callback = cb.affect!
51
52
setup = [
53
"GLM wave speed scaling" => glm_speed_callback.glm_scale,
54
"Expected CFL number" => glm_speed_callback.cfl,
55
"Selected semidiscretizations" => glm_speed_callback.semi_indices
56
]
57
summary_box(io, "GlmSpeedCallback", setup)
58
end
59
end
60
61
function GlmSpeedCallback(; glm_scale = 0.5, cfl, semi_indices = Int[])
62
@assert 0<=glm_scale<=1 "glm_scale must be between 0 and 1"
63
64
glm_speed_callback = GlmSpeedCallback{typeof(glm_scale), typeof(cfl)}(glm_scale,
65
cfl,
66
semi_indices)
67
68
DiscreteCallback(glm_speed_callback, glm_speed_callback, # the first one is the condition, the second the affect!
69
save_positions = (false, false),
70
initialize = initialize!)
71
end
72
73
function initialize!(cb::DiscreteCallback{Condition, Affect!}, u, t,
74
integrator) where {Condition, Affect! <: GlmSpeedCallback}
75
cb.affect!(integrator)
76
end
77
78
# this method is called to determine whether the callback should be activated
79
function (glm_speed_callback::GlmSpeedCallback)(u, t, integrator)
80
return true
81
end
82
83
function update_cleaning_speed!(semi, glm_speed_callback, dt, t)
84
@unpack glm_scale, cfl = glm_speed_callback
85
86
mesh, equations, solver, cache = mesh_equations_solver_cache(semi)
87
88
# compute time step for GLM linear advection equation with c_h=1 (redone due to the possible AMR)
89
if cfl isa Real # Case for constant CFL
90
c_h_deltat = calc_dt_for_cleaning_speed(cfl, mesh, equations, solver, cache)
91
else # Variable CFL
92
c_h_deltat = calc_dt_for_cleaning_speed(cfl(t), mesh, equations, solver, cache)
93
end
94
95
# c_h is proportional to its own time step divided by the complete MHD time step
96
# We use @reset here since the equations are immutable (to work on GPUs etc.).
97
# Thus, we need to modify the equations field of the semidiscretization.
98
@reset equations.c_h = glm_scale * c_h_deltat / dt
99
semi.equations = equations
100
101
return semi
102
end
103
104
# This method is called as callback after the StepsizeCallback during the time integration.
105
@inline function (glm_speed_callback::GlmSpeedCallback)(integrator)
106
dt = get_proposed_dt(integrator)
107
semi = integrator.p
108
109
# Call the appropriate update function (this indirection allows to specialize on,
110
# e.g., the semidiscretization type)
111
update_cleaning_speed!(semi, glm_speed_callback, dt, integrator.t)
112
113
# avoid re-evaluating possible FSAL stages
114
u_modified!(integrator, false)
115
116
return nothing
117
end
118
119
include("glm_speed_dg.jl")
120
end # @muladd
121
122