Path: blob/main/src/semidiscretization/semidiscretization_euler_gravity.jl
5586 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"""8ParametersEulerGravity(; background_density=0.0,9gravitational_constant=1.0,10cfl=1.0,11resid_tol=1.0e-4,12n_iterations_max=10^4,13timestep_gravity=timestep_gravity_erk52_3Sstar!)1415Set up parameters for the gravitational part of a [`SemidiscretizationEulerGravity`](@ref).1617# Arguments18- `background_density<:Real`: Constant background/reference density ρ₀ which is subtracted from the (Euler) density19in the RHS source term computation of the gravity solver.20- `gravitational_constant<:Real`: Gravitational constant G which needs to be in consistent units with the21density and velocity fields.22- `cfl<:Real`: CFL number used for the pseudo-time stepping to advance the hyperbolic diffusion equations into steady state.23- `resid_tol<:Real`: Absolute tolerance for the residual of the hyperbolic diffusion equations which are solved to24(approximately) steady state.25- `n_iterations_max::Int`: Maximum number of iterations of the pseudo-time gravity solver.26If `n_iterations <= 0` the solver will iterate until the residual is less or equal `resid_tol`.27This can cause an infinite loop if the solver does not converge!28- `timestep_gravity`: Function to advance the gravity solver by one pseudo-time step.29There are three optimized methods available:301) `timestep_gravity_erk51_3Sstar!` (first-order),312) `timestep_gravity_erk52_3Sstar!` (second-order),323) `timestep_gravity_erk53_3Sstar!` (third-order).33Additionally, `timestep_gravity_carpenter_kennedy_erk54_2N!` (fourth-order) can be used.34"""35struct ParametersEulerGravity{RealT <: Real, TimestepGravity}36background_density :: RealT # aka rho037gravitational_constant :: RealT # aka G38cfl :: RealT # CFL number for the gravity solver39resid_tol :: RealT # Hyp.-Diff. Eq. steady state tolerance40n_iterations_max :: Int # Max. number of iterations of the pseudo-time gravity solver41timestep_gravity :: TimestepGravity42end4344function ParametersEulerGravity(; background_density = 0.0,45gravitational_constant = 1.0,46cfl = 1.0,47resid_tol = 1.0e-4,48n_iterations_max = 10^4,49timestep_gravity = timestep_gravity_erk52_3Sstar!)50background_density, gravitational_constant, cfl, resid_tol = promote(background_density,51gravitational_constant,52cfl, resid_tol)53return ParametersEulerGravity(background_density, gravitational_constant, cfl,54resid_tol,55n_iterations_max, timestep_gravity)56end5758function Base.show(io::IO, parameters::ParametersEulerGravity)59@nospecialize parameters # reduce precompilation time6061print(io, "ParametersEulerGravity(")62print(io, "background_density=", parameters.background_density)63print(io, ", gravitational_constant=", parameters.gravitational_constant)64print(io, ", cfl=", parameters.cfl)65print(io, ", n_iterations_max=", parameters.n_iterations_max)66print(io, ", timestep_gravity=", parameters.timestep_gravity)67print(io, ")")68return nothing69end70function Base.show(io::IO, ::MIME"text/plain", parameters::ParametersEulerGravity)71@nospecialize parameters # reduce precompilation time7273if get(io, :compact, false)74show(io, parameters)75else76setup = [77"background density (ρ₀)" => parameters.background_density,78"gravitational constant (G)" => parameters.gravitational_constant,79"CFL (gravity)" => parameters.cfl,80"max. #iterations" => parameters.n_iterations_max,81"time integrator" => parameters.timestep_gravity82]83summary_box(io, "ParametersEulerGravity", setup)84end85end8687"""88SemidiscretizationEulerGravity8990A struct containing everything needed to describe a spatial semidiscretization91of a the compressible Euler equations with self-gravity, reformulating the92Poisson equation for the gravitational potential as steady-state problem of93the hyperbolic diffusion equations.94- Michael Schlottke-Lakemper, Andrew R. Winters, Hendrik Ranocha, Gregor J. Gassner (2020)95"A purely hyperbolic discontinuous Galerkin approach for self-gravitating gas dynamics"96[arXiv: 2008.10593](https://arXiv.org/abs/2008.10593)97"""98struct SemidiscretizationEulerGravity{SemiEuler, SemiGravity,99Parameters <: ParametersEulerGravity, Cache} <:100AbstractSemidiscretization101semi_euler :: SemiEuler102semi_gravity :: SemiGravity103parameters :: Parameters104cache :: Cache105performance_counter :: PerformanceCounter106gravity_counter :: PerformanceCounter107end108# We assume some properties of the fields of the semidiscretization, e.g.,109# the `equations` and the `mesh` should have the same dimension. We check these110# properties in the outer constructor defined below. While we could ensure111# them even better in an inner constructor, we do not use this approach to112# simplify the integration with Adapt.jl for GPU usage, see113# https://github.com/trixi-framework/Trixi.jl/pull/2677#issuecomment-3591789921114115"""116SemidiscretizationEulerGravity(semi_euler::SemiEuler, semi_gravity::SemiGravity, parameters)117118Construct a semidiscretization of the compressible Euler equations with self-gravity.119`parameters` should be given as [`ParametersEulerGravity`](@ref).120"""121function SemidiscretizationEulerGravity(semi_euler::SemiEuler,122semi_gravity::SemiGravity,123parameters) where124{Mesh,125SemiEuler <:126SemidiscretizationHyperbolic{Mesh, <:AbstractCompressibleEulerEquations},127SemiGravity <:128SemidiscretizationHyperbolic{Mesh, <:AbstractHyperbolicDiffusionEquations}}129@assert ndims(semi_euler) == ndims(semi_gravity)130@assert typeof(semi_euler.mesh) == typeof(semi_gravity.mesh)131@assert polydeg(semi_euler.solver) == polydeg(semi_gravity.solver)132133u_ode = compute_coefficients(zero(real(semi_gravity)), semi_gravity)134du_ode = similar(u_ode)135# Registers for gravity solver, tailored to the 2N and 3S* methods implemented below136u_tmp1_ode = similar(u_ode)137u_tmp2_ode = similar(u_ode)138cache = (; u_ode, du_ode, u_tmp1_ode, u_tmp2_ode)139140performance_counter = PerformanceCounter()141gravity_counter = PerformanceCounter()142143return SemidiscretizationEulerGravity{typeof(semi_euler), typeof(semi_gravity),144typeof(parameters), typeof(cache)}(semi_euler,145semi_gravity,146parameters,147cache,148performance_counter,149gravity_counter)150end151152function remake(semi::SemidiscretizationEulerGravity;153uEltype = real(semi.semi_gravity.solver),154semi_euler = semi.semi_euler,155semi_gravity = semi.semi_gravity,156parameters = semi.parameters)157semi_euler = remake(semi_euler, uEltype = uEltype)158semi_gravity = remake(semi_gravity, uEltype = uEltype)159160# Recreate cache, i.e., registers for u with e.g. AD datatype161u_ode = compute_coefficients(zero(real(semi_gravity)), semi_gravity)162du_ode = similar(u_ode)163u_tmp1_ode = similar(u_ode)164u_tmp2_ode = similar(u_ode)165cache = (; u_ode, du_ode, u_tmp1_ode, u_tmp2_ode)166167return SemidiscretizationEulerGravity(semi_euler, semi_gravity, parameters)168end169170function Base.show(io::IO, semi::SemidiscretizationEulerGravity)171@nospecialize semi # reduce precompilation time172173print(io, "SemidiscretizationEulerGravity using")174print(io, semi.semi_euler)175print(io, ", ", semi.semi_gravity)176print(io, ", ", semi.parameters)177print(io, ", cache(")178for (idx, key) in enumerate(keys(semi.cache))179idx > 1 && print(io, " ")180print(io, key)181end182print(io, "))")183return nothing184end185186function Base.show(io::IO, mime::MIME"text/plain", semi::SemidiscretizationEulerGravity)187@nospecialize semi # reduce precompilation time188189if get(io, :compact, false)190show(io, semi)191else192summary_header(io, "SemidiscretizationEulerGravity")193summary_line(io, "semidiscretization Euler",194semi.semi_euler |> typeof |> nameof)195show(increment_indent(io), mime, semi.semi_euler)196summary_line(io, "semidiscretization gravity",197semi.semi_gravity |> typeof |> nameof)198show(increment_indent(io), mime, semi.semi_gravity)199summary_line(io, "parameters", semi.parameters |> typeof |> nameof)200show(increment_indent(io), mime, semi.parameters)201summary_footer(io)202end203end204205# The compressible Euler semidiscretization is considered to be the main semidiscretization.206# The hyperbolic diffusion equations part is only used internally to update the gravitational207# potential during an rhs! evaluation of the flow solver.208@inline function mesh_equations_solver_cache(semi::SemidiscretizationEulerGravity)209return mesh_equations_solver_cache(semi.semi_euler)210end211212@inline Base.ndims(semi::SemidiscretizationEulerGravity) = ndims(semi.semi_euler)213214@inline Base.real(semi::SemidiscretizationEulerGravity) = real(semi.semi_euler)215216# computes the coefficients of the initial condition217@inline function compute_coefficients(t, semi::SemidiscretizationEulerGravity)218compute_coefficients!(semi.cache.u_ode, t, semi.semi_gravity)219return compute_coefficients(t, semi.semi_euler)220end221222# computes the coefficients of the initial condition and stores the Euler part in `u_ode`223@inline function compute_coefficients!(u_ode, t, semi::SemidiscretizationEulerGravity)224compute_coefficients!(semi.cache.u_ode, t, semi.semi_gravity)225return compute_coefficients!(u_ode, t, semi.semi_euler)226end227228@inline function calc_error_norms(func, u, t, analyzer,229semi::SemidiscretizationEulerGravity, cache_analysis)230return calc_error_norms(func, u, t, analyzer, semi.semi_euler, cache_analysis)231end232233# Coupled Euler and gravity solver at each Runge-Kutta stage,234# corresponding to Algorithm 2 in Schlottke-Lakemper et al. (2020),235# https://dx.doi.org/10.1016/j.jcp.2021.110467236function rhs!(du_ode, u_ode, semi::SemidiscretizationEulerGravity, t)237@unpack semi_euler, semi_gravity, cache = semi238239u_euler = wrap_array(u_ode, semi_euler)240du_euler = wrap_array(du_ode, semi_euler)241u_gravity = wrap_array(cache.u_ode, semi_gravity)242n_elements = size(u_euler)[end]243244time_start = time_ns()245246# standard semidiscretization of the compressible Euler equations247@trixi_timeit timer() "Euler solver" rhs!(du_ode, u_ode, semi_euler, t)248249# compute gravitational potential and forces250@trixi_timeit timer() "gravity solver" update_gravity!(semi, u_ode)251252# add gravitational source source_terms to the Euler part253if ndims(semi_euler) == 1254@threaded for element in 1:n_elements255@views @. du_euler[2, .., element] -= u_euler[1, .., element] *256u_gravity[2, .., element]257@views @. du_euler[3, .., element] -= u_euler[2, .., element] *258u_gravity[2, .., element]259end260elseif ndims(semi_euler) == 2261@threaded for element in 1:n_elements262@views @. du_euler[2, .., element] -= u_euler[1, .., element] *263u_gravity[2, .., element]264@views @. du_euler[3, .., element] -= u_euler[1, .., element] *265u_gravity[3, .., element]266@views @. du_euler[4, .., element] -= (u_euler[2, .., element] *267u_gravity[2, .., element] +268u_euler[3, .., element] *269u_gravity[3, .., element])270end271elseif ndims(semi_euler) == 3272@threaded for element in 1:n_elements273@views @. du_euler[2, .., element] -= u_euler[1, .., element] *274u_gravity[2, .., element]275@views @. du_euler[3, .., element] -= u_euler[1, .., element] *276u_gravity[3, .., element]277@views @. du_euler[4, .., element] -= u_euler[1, .., element] *278u_gravity[4, .., element]279@views @. du_euler[5, .., element] -= (u_euler[2, .., element] *280u_gravity[2, .., element] +281u_euler[3, .., element] *282u_gravity[3, .., element] +283u_euler[4, .., element] *284u_gravity[4, .., element])285end286else287error("Number of dimensions $(ndims(semi_euler)) not supported.")288end289290runtime = time_ns() - time_start291put!(semi.performance_counter, runtime)292293return nothing294end295296# TODO: Taal refactor, add some callbacks or so within the gravity update to allow investigating/optimizing it297function update_gravity!(semi::SemidiscretizationEulerGravity, u_ode)298@unpack semi_euler, semi_gravity, parameters, gravity_counter, cache = semi299300u_euler = wrap_array(u_ode, semi_euler)301u_gravity = wrap_array(cache.u_ode, semi_gravity)302du_gravity = wrap_array(cache.du_ode, semi_gravity)303304# set up main loop305finalstep = false306@unpack n_iterations_max, cfl, resid_tol, timestep_gravity = parameters307iter = 0308tau = zero(real(semi_gravity.solver)) # Pseudo-time309310# iterate gravity solver until convergence or maximum number of iterations are reached311@unpack equations = semi_gravity312while !finalstep313dtau = @trixi_timeit timer() "calculate dtau" begin314cfl * max_dt(u_gravity, tau, semi_gravity.mesh,315have_constant_speed(equations), equations,316semi_gravity.solver, semi_gravity.cache)317end318319# evolve solution by one pseudo-time step320time_start = time_ns()321timestep_gravity(cache, u_euler, tau, dtau, parameters, semi_gravity)322runtime = time_ns() - time_start323put!(gravity_counter, runtime)324325# update iteration counter326iter += 1327tau += dtau328329# check if we reached the maximum number of iterations330if n_iterations_max > 0 && iter >= n_iterations_max331@warn "Max iterations reached: Gravity solver failed to converge!" residual=maximum(abs,332@views du_gravity[1,333..,334:]) tau=tau dtau=dtau335finalstep = true336end337338# this is an absolute tolerance check339if maximum(abs, @views du_gravity[1, .., :]) <= resid_tol340finalstep = true341end342end343344return nothing345end346347# Integrate gravity solver for 2N-type low-storage schemes348function timestep_gravity_2N!(cache, u_euler, tau, dtau, gravity_parameters,349semi_gravity,350a, b, c)351G = gravity_parameters.gravitational_constant352rho0 = gravity_parameters.background_density353grav_scale = -4.0 * pi * G354355# Note that `u_ode` is `u_gravity` in `rhs!` above356@unpack u_ode, du_ode, u_tmp1_ode = cache357n_elements = size(u_euler)[end]358359u_tmp1_ode .= zero(eltype(u_tmp1_ode))360du_gravity = wrap_array(du_ode, semi_gravity)361362for stage in eachindex(c)363tau_stage = tau + dtau * c[stage]364365# rhs! has the source term for the harmonic problem366# We don't need a `@trixi_timeit timer() "rhs!"` here since that's already367# included in the `rhs!` call.368rhs!(du_ode, u_ode, semi_gravity, tau_stage)369370# Source term: Jeans instability OR coupling convergence test OR blast wave371# put in gravity source term proportional to Euler density372# OBS! subtract off the background density ρ_0 (spatial mean value)373# Note: Adding to `du_gravity` is essentially adding to `du_ode`!374@threaded for element in 1:n_elements375@views @. du_gravity[1, .., element] += grav_scale *376(u_euler[1, .., element] - rho0)377end378379a_stage = a[stage]380b_stage_dtau = b[stage] * dtau381@trixi_timeit timer() "Runge-Kutta step" begin382@threaded for idx in eachindex(u_ode)383u_tmp1_ode[idx] = du_ode[idx] - u_tmp1_ode[idx] * a_stage384u_ode[idx] += u_tmp1_ode[idx] * b_stage_dtau385end386end387end388389return nothing390end391392function timestep_gravity_carpenter_kennedy_erk54_2N!(cache, u_euler, tau, dtau,393gravity_parameters, semi_gravity)394# Coefficients for Carpenter's 5-stage 4th-order low-storage Runge-Kutta method395a = SVector(0.0,396567301805773.0 / 1357537059087.0,3972404267990393.0 / 2016746695238.0,3983550918686646.0 / 2091501179385.0,3991275806237668.0 / 842570457699.0)400b = SVector(1432997174477.0 / 9575080441755.0,4015161836677717.0 / 13612068292357.0,4021720146321549.0 / 2090206949498.0,4033134564353537.0 / 4481467310338.0,4042277821191437.0 / 14882151754819.0)405c = SVector(0.0,4061432997174477.0 / 9575080441755.0,4072526269341429.0 / 6820363962896.0,4082006345519317.0 / 3224310063776.0,4092802321613138.0 / 2924317926251.0)410411return timestep_gravity_2N!(cache, u_euler, tau, dtau, gravity_parameters,412semi_gravity, a, b, c)413end414415# Integrate gravity solver for 3S*-type low-storage schemes416function timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,417semi_gravity,418gamma1, gamma2, gamma3, beta, delta, c)419G = gravity_parameters.gravitational_constant420rho0 = gravity_parameters.background_density421grav_scale = -4 * G * pi422423# Note that `u_ode` is `u_gravity` in `rhs!` above424@unpack u_ode, du_ode, u_tmp1_ode, u_tmp2_ode = cache425n_elements = size(u_euler)[end]426427u_tmp1_ode .= zero(eltype(u_tmp1_ode))428u_tmp2_ode .= u_ode429du_gravity = wrap_array(du_ode, semi_gravity)430431for stage in eachindex(c)432tau_stage = tau + dtau * c[stage]433434# rhs! has the source term for the harmonic problem435# We don't need a `@trixi_timeit timer() "rhs!"` here since that's already436# included in the `rhs!` call.437rhs!(du_ode, u_ode, semi_gravity, tau_stage)438439# Source term: Jeans instability OR coupling convergence test OR blast wave440# put in gravity source term proportional to Euler density441# OBS! subtract off the background density ρ_0 around which the Jeans instability is perturbed442# Note: Adding to `du_gravity` is essentially adding to `du_ode`!443@threaded for element in 1:n_elements444@views @. du_gravity[1, .., element] += grav_scale *445(u_euler[1, .., element] - rho0)446end447448delta_stage = delta[stage]449gamma1_stage = gamma1[stage]450gamma2_stage = gamma2[stage]451gamma3_stage = gamma3[stage]452beta_stage_dtau = beta[stage] * dtau453@trixi_timeit timer() "Runge-Kutta step" begin454@threaded for idx in eachindex(u_ode)455# See Algorithm 1 (3S* method) in Schlottke-Lakemper et al. (2020)456u_tmp1_ode[idx] += delta_stage * u_ode[idx]457u_ode[idx] = (gamma1_stage * u_ode[idx] +458gamma2_stage * u_tmp1_ode[idx] +459gamma3_stage * u_tmp2_ode[idx] +460beta_stage_dtau * du_ode[idx])461end462end463end464465return nothing466end467468# First-order, 5-stage, 3S*-storage optimized method469function timestep_gravity_erk51_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,470semi_gravity)471# New 3Sstar coefficients optimized for polynomials of degree polydeg=3472# and examples/parameters_hypdiff_lax_friedrichs.toml473# 5 stages, order 1474gamma1 = SVector(0.0000000000000000E+00, 5.2910412316555866E-01,4752.8433964362349406E-01, -1.4467571130907027E+00,4767.5592215948661057E-02)477gamma2 = SVector(1.0000000000000000E+00, 2.6366970460864109E-01,4783.7423646095836322E-01, 7.8786901832431289E-01,4793.7754129043053775E-01)480gamma3 = SVector(0.0000000000000000E+00, 0.0000000000000000E+00,4810.0000000000000000E+00, 8.0043329115077388E-01,4821.3550099149374278E-01)483beta = SVector(1.9189497208340553E-01, 5.4506406707700059E-02,4841.2103893164085415E-01, 6.8582252490550921E-01,4858.7914657211972225E-01)486delta = SVector(1.0000000000000000E+00, 7.8593091509463076E-01,4871.2639038717454840E-01, 1.7726945920209813E-01,4880.0000000000000000E+00)489c = SVector(0.0000000000000000E+00, 1.9189497208340553E-01, 1.9580448818599061E-01,4902.4241635859769023E-01, 5.0728347557552977E-01)491492return timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,493semi_gravity,494gamma1, gamma2, gamma3, beta, delta, c)495end496497# Second-order, 5-stage, 3S*-storage optimized method498function timestep_gravity_erk52_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,499semi_gravity)500# New 3Sstar coefficients optimized for polynomials of degree polydeg=3501# and examples/parameters_hypdiff_lax_friedrichs.toml502# 5 stages, order 2503gamma1 = SVector(0.0000000000000000E+00, 5.2656474556752575E-01,5041.0385212774098265E+00, 3.6859755007388034E-01,505-6.3350615190506088E-01)506gamma2 = SVector(1.0000000000000000E+00, 4.1892580153419307E-01,507-2.7595818152587825E-02, 9.1271323651988631E-02,5086.8495995159465062E-01)509gamma3 = SVector(0.0000000000000000E+00, 0.0000000000000000E+00,5100.0000000000000000E+00, 4.1301005663300466E-01,511-5.4537881202277507E-03)512beta = SVector(4.5158640252832094E-01, 7.5974836561844006E-01,5133.7561630338850771E-01, 2.9356700007428856E-02,5142.5205285143494666E-01)515delta = SVector(1.0000000000000000E+00, 1.3011720142005145E-01,5162.6579275844515687E-01, 9.9687218193685878E-01,5170.0000000000000000E+00)518c = SVector(0.0000000000000000E+00, 4.5158640252832094E-01, 1.0221535725056414E+00,5191.4280257701954349E+00, 7.1581334196229851E-01)520521return timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,522semi_gravity,523gamma1, gamma2, gamma3, beta, delta, c)524end525526# Third-order, 5-stage, 3S*-storage optimized method527function timestep_gravity_erk53_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,528semi_gravity)529# New 3Sstar coefficients optimized for polynomials of degree polydeg=3530# and examples/parameters_hypdiff_lax_friedrichs.toml531# 5 stages, order 3532gamma1 = SVector(0.0000000000000000E+00, 6.9362208054011210E-01,5339.1364483229179472E-01, 1.3129305757628569E+00,534-1.4615811339132949E+00)535gamma2 = SVector(1.0000000000000000E+00, 1.3224582239681788E+00,5362.4213162353103135E-01, -3.8532017293685838E-01,5371.5603355704723714E+00)538gamma3 = SVector(0.0000000000000000E+00, 0.0000000000000000E+00,5390.0000000000000000E+00, 3.8306787039991996E-01,540-3.5683121201711010E-01)541beta = SVector(8.4476964977404881E-02, 3.0834660698015803E-01,5423.2131664733089232E-01, 2.8783574345390539E-01,5438.2199204703236073E-01)544delta = SVector(1.0000000000000000E+00, -7.6832695815481578E-01,5451.2497251501714818E-01, 1.4496404749796306E+00,5460.0000000000000000E+00)547c = SVector(0.0000000000000000E+00, 8.4476964977404881E-02, 2.8110631488732202E-01,5485.7093842145029405E-01, 7.2999896418559662E-01)549550return timestep_gravity_3Sstar!(cache, u_euler, tau, dtau, gravity_parameters,551semi_gravity,552gamma1, gamma2, gamma3, beta, delta, c)553end554555# TODO: Taal decide, where should specific parts like these be?556@inline function save_solution_file(u_ode, t, dt, iter,557semi::SemidiscretizationEulerGravity,558solution_callback,559element_variables = Dict{Symbol, Any}();560system = "")561# If this is called already as part of a multi-system setup (i.e., system is non-empty),562# we build a combined system name563if !isempty(system)564system_euler = system * "_euler"565system_gravity = system * "_gravity"566else567system_euler = "euler"568system_gravity = "gravity"569end570571u_euler = wrap_array_native(u_ode, semi.semi_euler)572filename_euler = save_solution_file(u_euler, t, dt, iter,573mesh_equations_solver_cache(semi.semi_euler)...,574solution_callback, element_variables,575system = system_euler)576577u_gravity = wrap_array_native(semi.cache.u_ode, semi.semi_gravity)578filename_gravity = save_solution_file(u_gravity, t, dt, iter,579mesh_equations_solver_cache(semi.semi_gravity)...,580solution_callback, element_variables,581system = system_gravity)582583return filename_euler, filename_gravity584end585586@inline function (amr_callback::AMRCallback)(u_ode,587semi::SemidiscretizationEulerGravity,588t, iter; kwargs...)589passive_args = ((semi.cache.u_ode,590mesh_equations_solver_cache(semi.semi_gravity)...),)591has_changed = amr_callback(u_ode, mesh_equations_solver_cache(semi.semi_euler)...,592semi, t, iter;593kwargs..., passive_args = passive_args)594595if has_changed596new_length = length(semi.cache.u_ode)597resize!(semi.cache.du_ode, new_length)598resize!(semi.cache.u_tmp1_ode, new_length)599resize!(semi.cache.u_tmp2_ode, new_length)600end601602return has_changed603end604end # @muladd605606607