Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix calc ch #353

Merged
merged 22 commits into from
Nov 23, 2020
Merged
Show file tree
Hide file tree
Changes from 16 commits
Commits
Show all changes
22 commits
Select commit Hold shift + click to select a range
a9b3caf
initial fix of the c_h computation a la Fluxo. Still ugly and involve…
andrewwinters5000 Nov 20, 2020
d509ef2
updated comments with known issue about c_h time sclaing. Tests will …
andrewwinters5000 Nov 20, 2020
114ba56
made glm speed update into a new callback. still needs cleanup and te…
andrewwinters5000 Nov 20, 2020
e440803
now callback is dimension agnostic
andrewwinters5000 Nov 20, 2020
291d413
glm callback now stores the cfl number for convenience. fixed bug whe…
andrewwinters5000 Nov 20, 2020
615fefa
removed commented code from stepsize callback
andrewwinters5000 Nov 20, 2020
1c6dd87
fixed bug in the order of arguments within the struct
andrewwinters5000 Nov 20, 2020
41a7562
added the glm callback to all the MHD examples
andrewwinters5000 Nov 20, 2020
97a71ee
all tests updated and pass locally
andrewwinters5000 Nov 20, 2020
5a7e631
address many comments and simplifications
andrewwinters5000 Nov 20, 2020
6cfb566
rename the new callback file
andrewwinters5000 Nov 20, 2020
13fe408
forgot to rename file in the include
andrewwinters5000 Nov 20, 2020
329e216
rename file
andrewwinters5000 Nov 20, 2020
4f7fe20
address comments, separate computation of ch for dg into a new file, …
andrewwinters5000 Nov 20, 2020
4d9b824
forgot to update how the timestep is extracted from the integrator
andrewwinters5000 Nov 20, 2020
4722d53
updated 3d tests to use new c_h calculation
andrewwinters5000 Nov 20, 2020
145afaf
slight changes to make glm callback structure match the other callbac…
andrewwinters5000 Nov 23, 2020
1edaabd
added assert error and more documentation for the GLM speed callback
andrewwinters5000 Nov 23, 2020
743866a
Merge branch 'master' into fix_calc_ch
andrewwinters5000 Nov 23, 2020
211cd60
added a check to ensure that the c_h is set to something in the ideal…
andrewwinters5000 Nov 23, 2020
e84ecd0
moved assert error to be before new memory is used
andrewwinters5000 Nov 23, 2020
01b8361
fixed error in assert and updated docstring
andrewwinters5000 Nov 23, 2020
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
15 changes: 9 additions & 6 deletions examples/2d/elixir_mhd_alfven_wave.jl
Original file line number Diff line number Diff line change
Expand Up @@ -36,22 +36,25 @@ analysis_callback = AnalysisCallback(semi, interval=analysis_interval, save_anal
extra_analysis_integrals=(entropy, energy_total,
energy_kinetic, energy_internal,
energy_magnetic, cross_helicity))

alive_callback = AliveCallback(analysis_interval=analysis_interval)

save_solution = SaveSolutionCallback(interval=10,
save_initial_solution=true,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.0)
cfl = 1.0
stepsize_callback = StepsizeCallback(cfl=cfl)

glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
stepsize_callback)

stepsize_callback,
glm_speed_callback)


###############################################################################
Expand Down
12 changes: 8 additions & 4 deletions examples/2d/elixir_mhd_alfven_wave_mortar.jl
Original file line number Diff line number Diff line change
Expand Up @@ -47,13 +47,17 @@ save_solution = SaveSolutionCallback(interval=10,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.0)
cfl = 1.0
stepsize_callback = StepsizeCallback(cfl=cfl)

glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
analysis_callback,
alive_callback,
save_solution,
stepsize_callback)
stepsize_callback,
glm_speed_callback)

###############################################################################
# run the simulation
Expand Down
14 changes: 9 additions & 5 deletions examples/2d/elixir_mhd_blast_wave.jl
Original file line number Diff line number Diff line change
Expand Up @@ -63,14 +63,18 @@ amr_callback = AMRCallback(semi, amr_controller,
adapt_initial_condition=true,
adapt_initial_condition_only_refine=true)

stepsize_callback = StepsizeCallback(cfl=0.3)
cfl = 0.3
stepsize_callback = StepsizeCallback(cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
amr_callback,
stepsize_callback)
amr_callback,
stepsize_callback,
glm_speed_callback)

###############################################################################
# run the simulation
Expand Down
12 changes: 8 additions & 4 deletions examples/2d/elixir_mhd_ec.jl
Original file line number Diff line number Diff line change
Expand Up @@ -41,13 +41,17 @@ save_solution = SaveSolutionCallback(interval=10,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.0)
cfl = 1.0
stepsize_callback = StepsizeCallback(cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
stepsize_callback)
stepsize_callback,
glm_speed_callback)


###############################################################################
Expand Down
14 changes: 9 additions & 5 deletions examples/2d/elixir_mhd_orszag_tang.jl
Original file line number Diff line number Diff line change
Expand Up @@ -63,14 +63,18 @@ amr_callback = AMRCallback(semi, amr_controller,
adapt_initial_condition=true,
adapt_initial_condition_only_refine=true)

stepsize_callback = StepsizeCallback(cfl=0.4)
cfl = 0.4
stepsize_callback = StepsizeCallback(cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
amr_callback,
stepsize_callback)
amr_callback,
stepsize_callback,
glm_speed_callback)

###############################################################################
# run the simulation
Expand Down
14 changes: 9 additions & 5 deletions examples/2d/elixir_mhd_rotor.jl
Original file line number Diff line number Diff line change
Expand Up @@ -64,14 +64,18 @@ amr_callback = AMRCallback(semi, amr_controller,
adapt_initial_condition=true,
adapt_initial_condition_only_refine=true)

stepsize_callback = StepsizeCallback(cfl=0.35)
cfl = 0.35
stepsize_callback = StepsizeCallback(cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
amr_callback,
stepsize_callback)
amr_callback,
stepsize_callback,
glm_speed_callback)

###############################################################################
# run the simulation
Expand Down
7 changes: 5 additions & 2 deletions examples/3d/elixir_mhd_alfven_wave.jl
Original file line number Diff line number Diff line change
Expand Up @@ -40,13 +40,16 @@ save_solution = SaveSolutionCallback(interval=10,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.4)
cfl = 1.4
stepsize_callback = StepsizeCallback(cfl=cfl)

glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback, alive_callback,
save_solution,
stepsize_callback)
stepsize_callback,
glm_speed_callback)


###############################################################################
Expand Down
8 changes: 6 additions & 2 deletions examples/3d/elixir_mhd_alfven_wave_mortar.jl
Original file line number Diff line number Diff line change
Expand Up @@ -43,12 +43,16 @@ save_solution = SaveSolutionCallback(interval=10,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.4)
cfl = 1.4
stepsize_callback = StepsizeCallback(cfl=cfl)

glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback, alive_callback,
save_solution,
stepsize_callback)
stepsize_callback,
glm_speed_callback)


###############################################################################
Expand Down
10 changes: 7 additions & 3 deletions examples/3d/elixir_mhd_ec.jl
Original file line number Diff line number Diff line change
Expand Up @@ -41,12 +41,16 @@ save_solution = SaveSolutionCallback(interval=10,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.4)
cfl = 1.4
stepsize_callback = StepsizeCallback(cfl=cfl)

glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

callbacks = CallbackSet(summary_callback,
analysis_callback, alive_callback,
analysis_callback, alive_callback,
save_solution,
stepsize_callback)
stepsize_callback,
glm_speed_callback)


###############################################################################
Expand Down
12 changes: 8 additions & 4 deletions examples/3d/elixir_mhd_orszag_tang.jl
Original file line number Diff line number Diff line change
Expand Up @@ -38,15 +38,19 @@ save_solution = SaveSolutionCallback(interval=10,
save_final_solution=true,
solution_variables=:primitive)

stepsize_callback = StepsizeCallback(cfl=1.1)
cfl = 1.1
stepsize_callback = StepsizeCallback(cfl=cfl)

glm_speed_callback = GlmSpeedCallback(glm_scale=0.5, cfl=cfl)

analysis_interval = 100

callbacks = CallbackSet(summary_callback,
analysis_callback,
callbacks = CallbackSet(summary_callback,
analysis_callback,
alive_callback,
save_solution,
stepsize_callback)
stepsize_callback,
glm_speed_callback)


###############################################################################
Expand Down
9 changes: 5 additions & 4 deletions src/Trixi.jl
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,8 @@ module Trixi
using LinearAlgebra: dot
using Printf: @printf, @sprintf, println

import DiffEqBase: ODEProblem, ODESolution, get_du, u_modified!, set_proposed_dt!, terminate!
import DiffEqBase: ODEProblem, ODESolution, get_du, u_modified!, set_proposed_dt!, terminate!,
get_proposed_dt
using DiffEqCallbacks: CallbackSet, DiscreteCallback
using EllipsisNotation # ..
using HDF5: h5open, attrs
Expand Down Expand Up @@ -107,9 +108,9 @@ export SemidiscretizationHyperbolic, semidiscretize, compute_coefficients, integ
export SemidiscretizationEulerGravity, ParametersEulerGravity,
timestep_gravity_erk52_3Sstar!, timestep_gravity_carpenter_kennedy_erk54_2N!

export SummaryCallback, SteadyStateCallback, AMRCallback, StepsizeCallback,
SaveRestartCallback, SaveSolutionCallback, AnalysisCallback, AliveCallback,
TrivialCallback
export SummaryCallback, SteadyStateCallback, AnalysisCallback, AliveCallback,
SaveRestartCallback, SaveSolutionCallback, AMRCallback, StepsizeCallback,
GlmSpeedCallback, TrivialCallback

export load_mesh, load_time

Expand Down
2 changes: 1 addition & 1 deletion src/auxiliary/precompile.jl
Original file line number Diff line number Diff line change
Expand Up @@ -230,6 +230,7 @@ function _precompile_manual_()
@assert Base.precompile(Tuple{Core.kwftype(typeof(Trixi.Type)),NamedTuple{(:analysis_interval,),Tuple{Int}},Type{AliveCallback}})
for RealT in (Float64,)
@assert Base.precompile(Tuple{Core.kwftype(typeof(Trixi.Type)),NamedTuple{(:cfl,),Tuple{RealT}},Type{StepsizeCallback}})
@assert Base.precompile(Tuple{Core.kwftype(typeof(Trixi.Type)),NamedTuple{(:glm_scale, :cfl),Tuple{RealT,RealT}},Type{GlmSpeedCallback}})
end
@assert Base.precompile(Tuple{Core.kwftype(typeof(Trixi.Type)),NamedTuple{(:interval, :save_final_restart),Tuple{Int,Bool}},Type{SaveRestartCallback}})
@assert Base.precompile(Tuple{Core.kwftype(typeof(Trixi.Type)),NamedTuple{(:interval, :save_initial_solution, :save_final_solution, :solution_variables),Tuple{Int,Bool,Bool,Symbol}},Type{SaveSolutionCallback}})
Expand Down Expand Up @@ -273,4 +274,3 @@ function _precompile_manual_()

return nothing
end

3 changes: 3 additions & 0 deletions src/callbacks_step/callbacks_step.jl
Original file line number Diff line number Diff line change
Expand Up @@ -36,6 +36,7 @@ end
# in a time step loop (however, callbacks are always executed *after* a step, thus it comes near
# the end here)
# * `StepsizeCallback` must come after AMR to accomodate potential changes in the minimum cell size
# * `GlmSpeedCallback` must come after computing time step size because it affects the value of c_h
include("summary.jl")
include("steady_state.jl")
include("analysis.jl")
Expand All @@ -44,6 +45,8 @@ include("save_restart.jl")
include("save_solution.jl")
include("amr.jl")
include("stepsize.jl")
include("glm_speed.jl")


# The `TrivialCallback` purposely does nothing: It allows to quickly disable specific callbacks
# when using `trixi_include` or `test_trixi_include`
Expand Down
70 changes: 70 additions & 0 deletions src/callbacks_step/glm_speed.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,70 @@

"""
GlmSpeedCallback(; glm_scale=0.5, cfl=1.0)

Update the divergence cleaning wave speed `c_h` according to the time step
computed in [`StepsizeCallback`](@ref) for the ideal GLM-MHD equations.
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
"""
struct GlmSpeedCallback{RealT<:Real}
glm_scale::RealT
cfl::RealT
end


function Base.show(io::IO, cb::DiscreteCallback{Condition,Affect!}) where {Condition, Affect!<:GlmSpeedCallback}
glm_speed_callback = cb.affect!
@unpack glm_scale, cfl = glm_speed_callback
print(io, "GlmSpeedCallback(glm_scale=", glm_scale, ", cfl=", cfl, ")")
end


function Base.show(io::IO, ::MIME"text/plain", cb::DiscreteCallback{Condition,Affect!}) where {Condition, Affect!<:GlmSpeedCallback}
if get(io, :compact, false)
show(io, cb)
else
glm_speed_callback = cb.affect!

setup = [
"GLM wave speed scaling" => glm_speed_callback.glm_scale,
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
"Expected CFL number" => glm_speed_callback.cfl,
]
summary_box(io, "GlmSpeedCallback", setup)
end
end


function GlmSpeedCallback(; glm_scale=0.5, cfl=1.0)
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
# when is the callback activated
condition = (u, t, integrator) -> true

glm_speed_callback = GlmSpeedCallback(glm_scale, cfl)

DiscreteCallback(condition, glm_speed_callback,
save_positions=(false,false),
initialize=initialize!)
end


function initialize!(cb::DiscreteCallback{Condition,Affect!}, u, t, integrator) where {Condition, Affect!<:GlmSpeedCallback}
cb.affect!(integrator)
end


# This method is called as callback after the StepsizeCallback during the time integration.
@inline function (glm_speed_callback::GlmSpeedCallback)(integrator)

dt = get_proposed_dt(integrator)
semi = integrator.p
_, equations, solver, cache = mesh_equations_solver_cache(semi)
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
@unpack glm_scale, cfl = glm_speed_callback

# compute time step for GLM linear advection equation with c_h=1 (redone due to the possible AMR)
c_h_deltat = calc_dt_for_cleaning_speed(cfl, equations, solver, cache)
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved

# c_h is proportional to its own time step divided by the complete MHD time step
equations.c_h = glm_scale * c_h_deltat / dt
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved

return nothing
end

include("glm_speed_dg.jl")
7 changes: 7 additions & 0 deletions src/callbacks_step/glm_speed_dg.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,7 @@

function calc_dt_for_cleaning_speed(cfl::Real, equations::AbstractIdealGlmMhdEquations, dg::DG, cache)
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
# compute time step for GLM linear advection equation with c_h=1 for the DG discretization on
# Cartesian meshes
max_scaled_speed_for_c_h = maximum(cache.elements.inverse_jacobian) * ndims(equations)
return cfl * 2 / (nnodes(dg) * max_scaled_speed_for_c_h)
andrewwinters5000 marked this conversation as resolved.
Show resolved Hide resolved
end
7 changes: 0 additions & 7 deletions src/callbacks_step/stepsize_dg2d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -5,13 +5,6 @@ function max_dt(u::AbstractArray{<:Any,4}, t, mesh::TreeMesh{2},
# e.g. for steady-state linear advection
max_scaled_speed = nextfloat(zero(t))

# FIXME: This should be implemented properly using another callback
# or something else, cf.
# https://github.com/trixi-framework/Trixi.jl/issues/257
if equations isa AbstractIdealGlmMhdEquations
equations.c_h = zero(equations.c_h)
end

for element in eachelement(dg, cache)
max_λ1 = max_λ2 = zero(max_scaled_speed)
for j in eachnode(dg), i in eachnode(dg)
Expand Down
8 changes: 0 additions & 8 deletions src/callbacks_step/stepsize_dg3d.jl
Original file line number Diff line number Diff line change
Expand Up @@ -5,13 +5,6 @@ function max_dt(u::AbstractArray{<:Any,5}, t, mesh::TreeMesh{3},
# e.g. for steady-state linear advection
max_scaled_speed = nextfloat(zero(t))

# FIXME: This should be implemented properly using another callback
# or something else, cf.
# https://github.com/trixi-framework/Trixi.jl/issues/257
if equations isa AbstractIdealGlmMhdEquations
equations.c_h = zero(equations.c_h)
end

for element in eachelement(dg, cache)
max_λ1 = max_λ2 = max_λ3 = zero(max_scaled_speed)
for k in eachnode(dg), j in eachnode(dg), i in eachnode(dg)
Expand Down Expand Up @@ -43,4 +36,3 @@ function max_dt(u::AbstractArray{<:Any,5}, t, mesh::TreeMesh{3},

return 2 / (nnodes(dg) * max_scaled_speed)
end

Loading