Evolution in Time

Once a Model is set up, it can be advanced in time with the evolve! function.

BoundaryLayerDynamics.evolve!Function
evolve!(model, tspan; dt, method, output)

Simulate a model’s evolution over time, optionally collecting data along the way.

Arguments

  • model::Model: The Model containing the current state as well as the discretized processes that describe the rate of change of the state.
  • tspan: The time span over which the evolution of the model should be simulated. Can be a single Number with the total duration or a Tuple of Numbers with the start and end times.

Keywords

  • dt: The (constant) size of the time step (required). Note that tspan has to be divisible by dt.
  • method = SSPRK33(): The time-integration method used to solve the semi-discretized initial-value problem as an ordinary differential equation.
  • output = []: A list of output modules that collect data during the simulation.
source

Time-Integration Methods

The time-integration methods are based on an extension of the standard formulation $\mathrm{d}\bm{\hat{s}}/\mathrm{d}t = f(\bm{\hat{s}}, t)$ that accommodates the role that pressure plays in incompressible flows. The methods solve a semi-discretized problem

\[\frac{\mathrm{d}\bm{\hat{s}}}{\mathrm{d}t} = f_r(\bm{\hat{s}}, t) + f_p(\bm{\hat{s}}, t)\]

where the right-hand-side term $f_r$ is computed directly from the state $\bm{\hat{s}}$ while $f_p$ is computed indirectly through a projection $P(\bm{\hat{s}})$ that can be used to enforce invariants of the state $\bm{\hat{s}}$.

Warning

None of the currently implemented processes are time-dependent, i.e. $f(\bm{\hat{s}}, t) = f(\bm{\hat{s}})$. It is likely that there are still errors in the handling of time-dependent processes and careful testing is highly recommended when implementing such processes.

Explicit Linear Multi-Step Methods

Explicit linear multi-step methods take the form

\[\bm{\hat{s}}^{(n+1)} = \sum_{i=1}^{m} \left( α_i \bm{\hat{s}}^{(n+1-i)} + Δt β_i f(\bm{\hat{s}}^{(n+1-i)}, t_i) \right) \,,\]

where each method is defined by its coefficients $m$, $α_i$, and $β_i$ and $\sum_{i=1}^m α_i = 1$ is required for consistency. The projection-based formulation of such a method is

\[\bm{\hat{s}}^{(n+1)} = P \left( \sum_{i=1}^{m} \left( α_i \bm{\hat{s}}^{(n+1-i)} + Δt β_i f_r(\bm{\hat{s}}^{(n+1-i)}, t_i) \right) \right) \,.\]

BoundaryLayerDynamics.ODEMethods.EulerType
Euler()

First-order explicit (“forward”) Euler method for time integration. This corresponds to an explicit linear multi-step methods with parameters $α₀=1$ and $β₀=1$.

source
BoundaryLayerDynamics.ODEMethods.AB2Type
AB2()

Second-order Adams–Bashforth method for time integration. This corresponds to an explicit linear multi-step methods with parameters $α₀=1$, $β₀=3/2$, and $β₁=−1/2$.

source

Explicit Runge–Kutta Methods

Explicit $m$-stage Runge–Kutta methods can be written as

\[\bm{\hat{s}}^{(n, i)} = \sum_{k=0}^{i-1} \left( α_{ik} \bm{\hat{s}}^{(n,k)} + Δt β_{ik} f\left(\bm{\hat{s}}^{(n,k)}, t^{(n,k)}\right) \right)\]

for $i=1,…,m$, with $\bm{\hat{s}}^{(n,0)} = \bm{\hat{s}}^{(n)}$ and $\bm{\hat{s}}^{(n+1)} = \bm{\hat{s}}^{(n,m)}$. $\sum_{k=0}^{i-1} α_{ik} = 1$ is required for consistency. The projection-based formulation of such a method is

\[\bm{\hat{s}}^{(n, i)} = P \left( \sum_{k=0}^{i-1} \left( α_{ik} \bm{\hat{s}}^{(n,k)} + Δt β_{ik} f_r\left(\bm{\hat{s}}^{(n,k)}, t^{(n,k)}\right) \right) \right) \,.\]

Output Modules

Output modules can collect information about a simulation while it is running. At each time step, they are provided with the state as well as any other fields that processes provide to the output modules by calling the log_sample! function.

BoundaryLayerDynamics.Logging.ProgressMonitorType
ProgressMonitor(; kwargs...)

Print information about an ongoing simulation to the standard output at regular intervals of wall time.

Keywords

  • frequency = 10: The wall-time interval (in seconds) between each progress report. Reports are always displayed at the end of a time step, if the wall time since the previous report exceeds the specified interval.
  • always_update = true: Print a dot after each time step between regular reports to indicate the ongoing progress of the simulation.
  • tstep = nothing: Passing the time step $Δt$ to the ProgressMonitor allows computing an advective Courant number $C = Δt / (∂₃^{max} |u|^{max})$. If not set, the advective time scale $1 / (∂₃^{max} |u|^{max})$ is reported instead.
source
BoundaryLayerDynamics.Logging.StepTimerType
StepTimer(; kwargs...)

Collect time stamps every $N$ steps to measure the compute time per time step and write it to a JSON file.

Keywords

  • path = "output/timestamps.json": The path to the output file.
  • frequency = 1: The number of time steps $N$ computed between each collected time stamp.
source
BoundaryLayerDynamics.Logging.MeanProfilesType
MeanProfiles(; kwargs...)

Collect profiles of terms, averaged over time and horizontal space. The time average is computed with the trapezoidal rule using samples collected before/after each time step. The averages are saved to HDF5 files, which can be written regularly during a simulation.

Keywords

  • profiles = (:vel1, :vel2, :vel3): Terms of which profiles are collected. Currently supported: :vel1, :vel2, :vel3, :adv1, :adv2, :adv3, :sgs11, :sgs12, :sgs13, :sgs22, :sgs23, :sgs33.
  • path = "output/profiles": Prefix of the path at which HDF5-files with the profiles are saved. The files are numbered, so the first file would be written to output/profiles-01.h5 for the default path.
  • output_frequency: The interval with which the output is written to disk. Note that this should be a multiple of the time step, otherwise some outputs may be skipped.
source
BoundaryLayerDynamics.Logging.SnapshotsType
Snapshots(; kwargs...)

Save the instantaneous state $\mathbf{\hat{s}}$ to a collection of Cartesian Binary Data files at regular intervals during the simulation.

Keywords

  • frequency: The frequency (in units of simulation time) at which snapshots should be saved. Note that this should be a multiple of the time step $Δt$.
  • path = "output/snapshots": Path to a folder inside which the snapshots will be saved.
  • centered = true: Save values at the $x₁$- and $x₂$-locations at the center of the intervals obtained by dividing the domain into $N₁×N₂$ segments. If set to false, the locations start at 0 instead and end one grid spacing before the end of the domain.
  • precision::DataType: The floating-point data type used in the output, either Float32 or Float64. By default this is set to the data type used for the simulation by the Model.
source