runge_kutta
moment_kinetics.runge_kutta
— ModuleRunge Kutta timestepping
moment_kinetics.runge_kutta.adaptive_timestep_update_t_params!
— Methodadaptive_timestep_update_t_params!(t_params, CFL_limits, error_norms,
total_points, error_norm_method, success,
nl_max_its_fraction, composition;
electron=false, local_max_dt::mk_float=Inf)
Use the calculated CFL_limits
and error_norms
to update the timestep in t_params
.
moment_kinetics.runge_kutta.local_error_norm
— Functionlocal_error_norm(error, f, rtol, atol)
local_error_norm(error, f, rtol, atol, neutral=false; method="Linf",
skip_r_inner=false, skip_z_lower=false, error_sum_zero=0.0)
Maximum error norm in the range owned by this MPI process, given by
\[\max(\frac{|\mathtt{error}|}{\mathtt{rtol}*|\mathtt{f}| + \mathtt{atol})\]
3 dimensional arrays (which represent moments) are treated as ion moments unless neutral=true
is passed.
method
can be "Linf" (to take the maximum error) or "L2" to take the root-mean-square (RMS) error.
skip_r_inner
and skip_z_lower
can be set to true to skip the contribution from the inner/lower boundaries, to avoid double-counting those points when using distributed-memory MPI.
error_sum_zero
should always have value 0.0, but is included so that different types can be used for L2sum. For testing, if we want consistency of results when using different numbers of processes (when the number of processes changes the order of operations in the sum is changed, which changes the rounding errors) then we have to use higher precision (i.e. use the Float128 type from the Quadmath package). The type of a 0.0 value can be set according to the high_precision_error_sum
option in the [timestepping]
section, and stored in a template-typed value in the t_params
object - when that value is passed in as the argument to error_sum_zero
, that type will be used for L2sum, and the type will be known at compile time, allowing this function to be efficient.
moment_kinetics.runge_kutta.rk_loworder_solution!
— MethodCalculate a lower-order approximation for the variable named var_symbol
, which can be used to calculate an error estimate for adaptive timestepping methods.
The lower-order approximation is stored in var_symbol
in scratch[2]
(as this entry should not be needed again after the lower-order approximation is calculated).
moment_kinetics.runge_kutta.rk_update_evolved_moments!
— Methoduse Runge Kutta to update any ion velocity moments evolved separately from the pdf
moment_kinetics.runge_kutta.rk_update_evolved_moments_electron!
— Methoduse Runge Kutta to update any electron velocity moments evolved separately from the pdf
moment_kinetics.runge_kutta.rk_update_evolved_moments_neutral!
— Methoduse Runge Kutta to update any neutral-particle velocity moments evolved separately from the pdf
moment_kinetics.runge_kutta.rk_update_variable!
— MethodUpdate the variable named var_symbol
in scratch
to the current Runge-Kutta stage istage
. The current value in scratch[istage+1]
is the result of the forward-Euler update, which needs to be corrected using values from previous stages with the Runge-Kutta coefficients. scratch_implicit
contains the results of backward-Euler updates, which are needed for IMEX timestepping schemes.
moment_kinetics.runge_kutta.setup_runge_kutta_coefficients!
— Methodgiven the number of Runge Kutta stages that are requested, returns the needed Runge Kutta coefficients; e.g., if f is the function to be updated, then f^{n+1}[stage+1] = rkcoef[1,stage]*f^{n} + rkcoef[2,stage]f^{n+1}[stage] + rk_coef[3,stage](f^{n}+dt*G[f^{n+1}[stage]]