Path: blob/main/src/time_integration/methods_SSP.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# Abstract base type for time integration schemes of explicit strong stability-preserving (SSP)8# Runge-Kutta (RK) methods. They are high-order time discretizations that guarantee the TVD property.9abstract type SimpleAlgorithmSSP <: AbstractTimeIntegrationAlgorithm end1011"""12SimpleSSPRK33(; stage_callbacks=())1314The third-order SSP Runge-Kutta method of Shu and Osher.1516## References1718- Shu, Osher (1988)19"Efficient Implementation of Essentially Non-oscillatory Shock-Capturing Schemes" (Eq. 2.18)20[DOI: 10.1016/0021-9991(88)90177-5](https://doi.org/10.1016/0021-9991(88)90177-5)21"""22struct SimpleSSPRK33{StageCallbacks} <: SimpleAlgorithmSSP23numerator_a::SVector{3, Float64}24numerator_b::SVector{3, Float64}25denominator::SVector{3, Float64}26c::SVector{3, Float64}27stage_callbacks::StageCallbacks2829function SimpleSSPRK33(; stage_callbacks = ())30# Mathematically speaking, it is not necessary for the algorithm to split the factors31# into numerator and denominator. Otherwise, however, rounding errors of the order of32# the machine accuracy will occur, which will add up over time and thus endanger the33# conservation of the simulation.34# See also https://github.com/trixi-framework/Trixi.jl/pull/1640.35numerator_a = SVector(0.0, 3.0, 1.0) # a = numerator_a / denominator36numerator_b = SVector(1.0, 1.0, 2.0) # b = numerator_b / denominator37denominator = SVector(1.0, 4.0, 3.0)38c = SVector(0.0, 1.0, 1 / 2)3940# Butcher tableau41# c | A42# 0 |43# 1 | 144# 1/2 | 1/4 1/445# --------------------46# b | 1/6 1/6 2/34748return new{typeof(stage_callbacks)}(numerator_a, numerator_b, denominator, c,49stage_callbacks)50end51end5253# This struct is needed to fake https://github.com/SciML/OrdinaryDiffEq.jl/blob/0c2048a502101647ac35faabd80da8a5645beac7/src/integrators/type.jl#L154mutable struct SimpleIntegratorSSPOptions{Callback, TStops}55callback::Callback # callbacks; used in Trixi56const adaptive::Bool # whether the algorithm is adaptive; ignored57dtmax::Float64 # ignored58const maxiters::Int # maximal number of time steps59tstops::TStops # tstops from https://diffeq.sciml.ai/v6.8/basics/common_solver_opts/#Output-Control-1; ignored60end6162function SimpleIntegratorSSPOptions(callback, tspan; maxiters = typemax(Int), kwargs...)63tstops_internal = BinaryHeap{eltype(tspan)}(FasterForward())64# We add last(tspan) to make sure that the time integration stops at the end time65push!(tstops_internal, last(tspan))66# We add 2 * last(tspan) because add_tstop!(integrator, t) is only called by DiffEqCallbacks.jl if tstops contains a time that is larger than t67# (https://github.com/SciML/DiffEqCallbacks.jl/blob/025dfe99029bd0f30a2e027582744528eb92cd24/src/iterative_and_periodic.jl#L92)68push!(tstops_internal, 2 * last(tspan))69return SimpleIntegratorSSPOptions{typeof(callback), typeof(tstops_internal)}(callback,70false,71Inf,72maxiters,73tstops_internal)74end7576# This struct is needed to fake https://github.com/SciML/OrdinaryDiffEq.jl/blob/0c2048a502101647ac35faabd80da8a5645beac7/src/integrators/type.jl#L7777# This implements the interface components described at78# https://diffeq.sciml.ai/v6.8/basics/integrator/#Handing-Integrators-179# which are used in Trixi.80mutable struct SimpleIntegratorSSP{RealT <: Real, uType,81Params, Sol, F, Alg,82SimpleIntegratorSSPOptions} <: AbstractTimeIntegrator83u::uType84du::uType85u_tmp::uType86t::RealT87tdir::RealT # DIRection of time integration, i.e., if one marches forward or backward in time88dt::RealT # current time step89dtcache::RealT # manually set time step90iter::Int # current number of time steps (iteration)91p::Params # will be the semidiscretization from Trixi92sol::Sol # faked93const f::F # `rhs!` of the semidiscretization94const alg::Alg # SimpleSSPRK3395opts::SimpleIntegratorSSPOptions96finalstep::Bool # added for convenience97const dtchangeable::Bool98const force_stepfail::Bool99end100101"""102add_tstop!(integrator::SimpleIntegratorSSP, t)103Add a time stop during the time integration process.104This function is called after the periodic SaveSolutionCallback to specify the next stop to save the solution.105"""106function add_tstop!(integrator::SimpleIntegratorSSP, t)107integrator.tdir * (t - integrator.t) < zero(integrator.t) &&108error("Tried to add a tstop that is behind the current time. This is strictly forbidden")109# We need to remove the first entry of tstops when a new entry is added.110# Otherwise, the simulation gets stuck at the previous tstop and dt is adjusted to zero.111if length(integrator.opts.tstops) > 1112pop!(integrator.opts.tstops)113end114return push!(integrator.opts.tstops, integrator.tdir * t)115end116117has_tstop(integrator::SimpleIntegratorSSP) = !isempty(integrator.opts.tstops)118first_tstop(integrator::SimpleIntegratorSSP) = first(integrator.opts.tstops)119120function init(ode::ODEProblem, alg::SimpleAlgorithmSSP;121dt, callback::Union{CallbackSet, Nothing} = nothing, kwargs...)122u = copy(ode.u0)123du = similar(u)124u_tmp = similar(u)125t = first(ode.tspan)126tdir = sign(ode.tspan[end] - ode.tspan[1])127iter = 0128integrator = SimpleIntegratorSSP(u, du, u_tmp, t, tdir, dt, dt, iter, ode.p,129(prob = ode,), ode.f, alg,130SimpleIntegratorSSPOptions(callback, ode.tspan;131kwargs...),132false, true, false)133134# Standard callbacks135initialize_callbacks!(callback, integrator)136137# Addition for `SimpleAlgorithmSSP` which may have stage callbacks138for stage_callback in alg.stage_callbacks139init_callback(stage_callback, integrator.p)140end141142return integrator143end144145function solve!(integrator::SimpleIntegratorSSP)146@unpack prob = integrator.sol147148integrator.finalstep = false149150@trixi_timeit timer() "main loop" while !integrator.finalstep151step!(integrator)152end153154# Empty the tstops array.155# This cannot be done in terminate!(integrator::SimpleIntegratorSSP) because DiffEqCallbacks.PeriodicCallbackAffect would return at error.156extract_all!(integrator.opts.tstops)157158for stage_callback in integrator.alg.stage_callbacks159finalize_callback(stage_callback, integrator.p)160end161162finalize_callbacks(integrator)163164return TimeIntegratorSolution((first(prob.tspan), integrator.t),165(prob.u0, integrator.u), prob)166end167168function step!(integrator::SimpleIntegratorSSP)169@unpack prob = integrator.sol170@unpack alg = integrator171t_end = last(prob.tspan)172callbacks = integrator.opts.callback173174@assert !integrator.finalstep175if isnan(integrator.dt)176error("time step size `dt` is NaN")177end178179modify_dt_for_tstops!(integrator)180181limit_dt!(integrator, t_end)182183@. integrator.u_tmp = integrator.u184for stage in eachindex(alg.c)185t_stage = integrator.t + integrator.dt * alg.c[stage]186# compute du187integrator.f(integrator.du, integrator.u, integrator.p, t_stage)188189# perform forward Euler step190@. integrator.u = integrator.u + integrator.dt * integrator.du191192for stage_callback in alg.stage_callbacks193stage_callback(integrator.u, integrator, stage)194end195196# perform convex combination197@. integrator.u = (alg.numerator_a[stage] * integrator.u_tmp +198alg.numerator_b[stage] * integrator.u) /199alg.denominator[stage]200end201integrator.iter += 1202integrator.t += integrator.dt203204@trixi_timeit timer() "Step-Callbacks" handle_callbacks!(callbacks, integrator)205206check_max_iter!(integrator)207208return nothing209end210211# get a cache where the RHS can be stored212get_tmp_cache(integrator::SimpleIntegratorSSP) = (integrator.u_tmp,)213214# some algorithms from DiffEq like FSAL-ones need to be informed when a callback has modified u215u_modified!(integrator::SimpleIntegratorSSP, ::Bool) = false216217# stop the time integration218function terminate!(integrator::SimpleIntegratorSSP)219integrator.finalstep = true220221return nothing222end223224"""225modify_dt_for_tstops!(integrator::SimpleIntegratorSSP)226Modify the time-step size to match the time stops specified in integrator.opts.tstops.227To avoid adding OrdinaryDiffEq to Trixi's dependencies, this routine is a copy of228https://github.com/SciML/OrdinaryDiffEq.jl/blob/d76335281c540ee5a6d1bd8bb634713e004f62ee/src/integrators/integrator_utils.jl#L38-L54229"""230function modify_dt_for_tstops!(integrator::SimpleIntegratorSSP)231if has_tstop(integrator)232tdir_t = integrator.tdir * integrator.t233tdir_tstop = first_tstop(integrator)234if integrator.opts.adaptive235integrator.dt = integrator.tdir *236min(abs(integrator.dt), abs(tdir_tstop - tdir_t)) # step! to the end237elseif iszero(integrator.dtcache) && integrator.dtchangeable238integrator.dt = integrator.tdir * abs(tdir_tstop - tdir_t)239elseif integrator.dtchangeable && !integrator.force_stepfail240# always try to step! with dtcache, but lower if a tstop241# however, if force_stepfail then don't set to dtcache, and no tstop worry242integrator.dt = integrator.tdir *243min(abs(integrator.dtcache), abs(tdir_tstop - tdir_t)) # step! to the end244end245end246247return nothing248end249250# used for AMR251function Base.resize!(integrator::SimpleIntegratorSSP, new_size)252resize!(integrator.u, new_size)253resize!(integrator.du, new_size)254resize!(integrator.u_tmp, new_size)255256return nothing257end258end # @muladd259260261