Path: blob/main/src/callbacks_step/glm_speed.jl
2055 views
# By default, Julia/LLVM does not use fused multiply-add operations (FMAs).1# Since these FMAs can increase the performance of many numerical algorithms,2# we need to opt-in explicitly.3# See https://ranocha.de/blog/Optimizing_EC_Trixi for further details.4@muladd begin5#! format: noindent67"""8GlmSpeedCallback(; glm_scale=0.5, cfl, semi_indices=())910Update the divergence cleaning wave speed `c_h` according to the time step11computed in [`StepsizeCallback`](@ref) for the ideal GLM-MHD equations, the multi-component12GLM-MHD equations, and the multi-ion GLM-MHD equations.13The `cfl` number should be set to the same value as for the time step size calculation.14As for standard [`StepsizeCallback`](@ref) `cfl` can be either set to a `Real` number or15a function of time `t` returning a `Real` number.1617The `glm_scale` ensures that the GLM wave speed is lower than the fastest physical waves in the MHD18solution and should thus be set to a value within the interval [0,1]. Note that `glm_scale = 0`19deactivates the divergence cleaning.2021In case of coupled semidiscretizations, specify for which `semi_index`, i.e. index of the22semidiscretization, the divergence cleaning should be applied. See also23[`SemidiscretizationCoupled`](@ref).24Note: `SemidiscretizationCoupled` and all related features are considered experimental and25may change at any time.26"""27struct GlmSpeedCallback{RealT <: Real, CflType}28glm_scale::RealT29cfl::CflType30semi_indices::Vector{Int}31end3233function Base.show(io::IO, cb::DiscreteCallback{<:Any, <:GlmSpeedCallback})34@nospecialize cb # reduce precompilation time3536glm_speed_callback = cb.affect!37@unpack glm_scale, cfl, semi_indices = glm_speed_callback38print(io, "GlmSpeedCallback(glm_scale=", glm_scale, ", cfl=", cfl, "semi_indices=",39semi_indices, ")")40end4142function Base.show(io::IO, ::MIME"text/plain",43cb::DiscreteCallback{<:Any, <:GlmSpeedCallback})44@nospecialize cb # reduce precompilation time4546if get(io, :compact, false)47show(io, cb)48else49glm_speed_callback = cb.affect!5051setup = [52"GLM wave speed scaling" => glm_speed_callback.glm_scale,53"Expected CFL number" => glm_speed_callback.cfl,54"Selected semidiscretizations" => glm_speed_callback.semi_indices55]56summary_box(io, "GlmSpeedCallback", setup)57end58end5960function GlmSpeedCallback(; glm_scale = 0.5, cfl, semi_indices = Int[])61@assert 0<=glm_scale<=1 "glm_scale must be between 0 and 1"6263glm_speed_callback = GlmSpeedCallback{typeof(glm_scale), typeof(cfl)}(glm_scale,64cfl,65semi_indices)6667DiscreteCallback(glm_speed_callback, glm_speed_callback, # the first one is the condition, the second the affect!68save_positions = (false, false),69initialize = initialize!)70end7172function initialize!(cb::DiscreteCallback{Condition, Affect!}, u, t,73integrator) where {Condition, Affect! <: GlmSpeedCallback}74cb.affect!(integrator)75end7677# this method is called to determine whether the callback should be activated78function (glm_speed_callback::GlmSpeedCallback)(u, t, integrator)79return true80end8182function update_cleaning_speed!(semi, glm_speed_callback, dt, t)83@unpack glm_scale, cfl = glm_speed_callback8485mesh, equations, solver, cache = mesh_equations_solver_cache(semi)8687# compute time step for GLM linear advection equation with c_h=1 (redone due to the possible AMR)88if cfl isa Real # Case for constant CFL89c_h_deltat = calc_dt_for_cleaning_speed(cfl, mesh, equations, solver, cache)90else # Variable CFL91c_h_deltat = calc_dt_for_cleaning_speed(cfl(t), mesh, equations, solver, cache)92end9394# c_h is proportional to its own time step divided by the complete MHD time step95# We use @reset here since the equations are immutable (to work on GPUs etc.).96# Thus, we need to modify the equations field of the semidiscretization.97@reset equations.c_h = glm_scale * c_h_deltat / dt98semi.equations = equations99100return semi101end102103# This method is called as callback after the StepsizeCallback during the time integration.104@inline function (glm_speed_callback::GlmSpeedCallback)(integrator)105dt = get_proposed_dt(integrator)106semi = integrator.p107108# Call the appropriate update function (this indirection allows to specialize on,109# e.g., the semidiscretization type)110update_cleaning_speed!(semi, glm_speed_callback, dt, integrator.t)111112# avoid re-evaluating possible FSAL stages113u_modified!(integrator, false)114115return nothing116end117118include("glm_speed_dg.jl")119end # @muladd120121122