FixedMeshRefinement

Documentation for FixedMeshRefinement.

FixedMeshRefinement.GridType
Grid

A struct representing the entire FMR grid, which consists of multiple Levels.

Fields

  • levels::Vector{Level}: A vector of Level objects.
  • base_dt::Float64: The time step of the coarsest level.
  • t::Float64: The current time of the simulation.
  • subcycling::Bool: A boolean indicating whether subcycling in time is enabled.
source
FixedMeshRefinement.LevelType
Level

A struct representing a single refinement level in the mesh.

Fields

  • index::Int: Index of the level.
  • num_state_variables::Int: Number of state variables.
  • num_diagnostic_variables::Int: Number of diagnostic variables.
  • num_tmp_variables::Int: Number of temporary variables.
  • num_interior_points::Int: Number of interior grid points.
  • num_ghost_points::Int: Number of ghost points on each side.
  • num_buffer_points::Int: Number of buffer points on each side for inter-level communication.
  • num_transition_points::NTuple{2,Int}: Number of points in the transition zone for mesh refinement.
  • num_boundary_points::NTuple{2,Int}: Number of boundary points on each side (ghost or buffer).
  • time_interpolation_order::Int: Order of time interpolation.
  • spatial_interpolation_order::Int: Order of spatial interpolation.
  • domain_box::Tuple{Float64,Float64}: The computational domain (interior) of this level.
  • physical_domain_box::Tuple{Float64,Float64}: The physical domain of the entire grid.
  • is_physical_boundary::NTuple{2,Bool}: Indicates if the level boundaries are physical boundaries.
  • dx::Float64: Grid spacing.
  • dt::Float64: Time step.
  • t::Float64: Current time of this level.
  • parent_indices::UnitRange{Int}: Range of indices of the parent level that this level covers.
  • offset_indices::UnitRange{Int}: Range of indices for OffsetArrays of this level.
source
Base.showMethod
Base.show(io::IO, grid::Grid)

Display a compact summary of the Grid.

source
Base.showMethod
Base.show(io::IO, level::Level)

Display a compact summary of the Level.

source
Base.showMethod
Base.show(io::IO, ::MIME"text/plain", grid::Grid)

Display a detailed summary of the Grid structure.

source
Base.showMethod
Base.show(io::IO, ::MIME"text/plain", level::Level)

Display a detailed summary of the Level.

source
FixedMeshRefinement.apply_transition_zone!Method
apply_transition_zone!(grid::Grid, l::Int, interp_in_time::Bool)

Apply a transition zone to smoothly connect coarse and fine grid solutions at the boundaries of a refinement level. This function blends the fine grid solution with a prolonged solution from the coarse grid.

Arguments

  • grid::Grid: The grid structure.
  • l::Int: The fine level index.
  • interp_in_time::Bool: Whether to use time interpolation.
source
FixedMeshRefinement.calc_Yn_from_kcs!Method
calc_Yn_from_kcs!(Yn_buffer, yn, kcs, dtc, interp_in_time, dytmp)

Calculate the intermediate Y values needed for Mongwane's subcycling method. These values represent the state on the coarse grid at intermediate times required by the fine grid's time stepping.

Arguments

  • Yn_buffer: Buffer to store the calculated Y values.
  • yn: State on the coarse grid at the beginning of the step.
  • kcs: The RK4 stages from the coarse grid time step.
  • dtc::Float64: The coarse grid time step.
  • interp_in_time::Bool: Flag indicating if interpolation in time is needed.
  • dytmp: Temporary storage for derivatives.
source
FixedMeshRefinement.cycle_state!Method
cycle_state!(level::Level)

Shift the time levels of the state variables in a Level. state[i] becomes state[i+1]. This is used to advance the solution in time.

source
FixedMeshRefinement.durationstringMethod
durationstring(nsec)

Convert a duration in seconds to a human-readable string format.

Arguments

  • nsec: Number of seconds to convert

Returns

  • A string representing the duration in the format:
    • "HH:MM:SS" for durations less than a day
    • "N days, HH:MM:SS" for durations between 1-9 days
    • "X.XX days" for durations greater than 9 days
source
FixedMeshRefinement.fill_buffer!Method
fill_buffer!(u, level::Level, stage::Int)

Fill the ghost cells of the solution array u from the Yn_buffer at a specific Runge-Kutta stage. This is part of the subcycling-in-time method. After copying, the buffer is filled with NaN to avoid accidental reuse.

Arguments

  • u: The solution array with ghost cells to be filled.
  • level::Level: The grid level.
  • stage::Int: The RK stage index.
  • fill_with_nan::Bool: If true, fill the buffer with NaN to avoid accidental reuse. Defaults to true.
source
FixedMeshRefinement.fine_to_coarse_indexMethod
fine_to_coarse_index(fine_level::Level, fidx::Int) -> Int

Convert a fine grid index fidx to the corresponding coarse grid index. This function is not defined for the base level. An error is thrown if the fine grid point does not align with any coarse grid point.

source
FixedMeshRefinement.get_stateFunction
get_state(level::Level, i::Int=0)

Return the state variables of the level as an OffsetArray. The optional argument i specifies the time level, where i=0 corresponds to the current time level n, i=-1 to n-1, etc.

source
FixedMeshRefinement.integrate_simpsonMethod
integrate_simpson(grid::Grid, getter::Function)

Integrate a variable over all levels of the grid using Simpson's rule. The getter function specifies which data to retrieve the data of the variable from each level. getter(level) should return an SubArray of OffsetArray.

source
FixedMeshRefinement.integrate_trapezoidalMethod
integrate_trapezoidal(grid::Grid, getter::Function)

Integrate a variable over all levels of the grid using the trapezoidal rule. The getter function specifies which data to retrieve the data of the variable from each level. getter(level) should return an SubArray of OffsetArray.

source
FixedMeshRefinement.merge_grid_levelsMethod
merge_grid_levels(grid::Grid, getter::Function)

Merge signle variable from all levels of the grid into a BlockArray. The getter function specifies which data to retrieve the data of the variable from each level. getter(level) should return a SubArray of OffsetArray from which data is read.

source
FixedMeshRefinement.prolongate!Method
prolongate!(grid::Grid, l::Int, interp_in_time::Bool)

Perform prolongation from a coarse grid (l-1) to a fine grid (l). This function fills the ghost cells of the fine grid by interpolating the data from the coarse grid in space and, optionally, in time.

Arguments

  • grid::Grid: The grid structure.
  • l::Int: The fine level index.
  • interp_in_time::Bool: Whether to use time interpolation.
source
FixedMeshRefinement.prolongate_mongwane!Method
prolongate_mongwane!(grid::Grid, l::Int, interp_in_time::Bool)

Perform prolongation from a coarse grid (l-1) to a fine grid (l) using Mongwane's subcycling-in-time method. This is used to set the ghost cell data for the fine grid when subcycling is enabled.

Arguments

  • grid::Grid: The grid structure.
  • l::Int: The fine level index.
  • interp_in_time::Bool: Whether to use time interpolation for the coarse grid state.
source
FixedMeshRefinement.restrict_injection!Method
restrict_injection!(grid::Grid, l::Int; apply_trans_zone=false)

Perform restriction from a fine grid (l+1) to a coarse grid (l). This operation transfers data from a finer grid to a coarser grid. Currently, it uses injection by taking the values from the fine grid at the locations corresponding to the coarse grid points.

Arguments

  • grid::Grid: The grid structure.
  • l::Int: The coarse level index.
  • apply_trans_zone::Bool: If true, the restriction is only applied outside of the transition zones. Defaults to false.
source
FixedMeshRefinement.rk4!Method
rk4!(level::Level, f::Function, p; mongwane::Bool=false)

Perform a single step of the classic 4th-order Runge-Kutta (RK4) method on a given level. It advances the solution by one time step dt.

Arguments

  • level::Level: The grid level to be updated.
  • f::Function: The function that computes the right-hand side of the ODEs. It should have the signature f(level, k, u, p, t).
  • p: Parameters to be passed to the RHS function f.
  • mongwane::Bool: If true, enables special buffer filling for Mongwane's subcycling method. Defaults to false.
source
FixedMeshRefinement.rk4_dense_output_y!Method
rk4_dense_output_y!(y, theta, h, yn, k)

Calculate the dense output for a 4th-order Runge-Kutta method at a fractional time theta within a time step h. This is used for interpolation in time within a single time step.

Arguments

  • y: The output array for the interpolated state.
  • theta::Float64: The fractional time within the interval [0, 1].
  • h::Float64: The time step size.
  • yn: The state at the beginning of the time step.
  • k: A vector of RK4 stages k1 through k4.
source
FixedMeshRefinement.shift_grid_boundaries!Method
shift_grid_boundaries!(grid::Grid, num_shift_points::NTuple{2,Int}; func_fill_extended::Function=fill_extended_grid_extrapolate!)

Shift the boundary of the grid by a given number of points; as a prerequisite, all levels must align with the physical boundary. If the number of shift points is positive, the grid will be extended, otherwise it will be shrunk. The tuple num_shift_points is the number of points to shift on the left and right boundaries, respectively.

Arguments

  • grid::Grid: The grid to shift boundaries of
  • num_shift_points::NTuple{2,Int}: Number of points to shift on left and right boundaries
  • func_fill_extended::Function=fill_extended_grid_extrapolate!: Function to compute fill values, takes index and level as arguments
source
FixedMeshRefinement.shift_level_boundaries!Method
shift_level_boundaries!(level::Level, num_shift_points::NTuple{2,Int}; func_fill_extended::Function=fill_extended_grid_extrapolate!)

Shift the boundary of a Level by a given number of points, as a prerequisite, the level must align with the physical boundary.

source
FixedMeshRefinement.spatial_interpolate!Method
spatial_interpolate!(res, u, i, order)

Performs spatial interpolation for prolongation using a set of predefined stencils of a given order. The result is stored in res. The interpolation uses values from u around index i.

Arguments

  • res: The destination array for the interpolated values.
  • u: The source data array.
  • i::Int: The index in u to center the interpolation stencil.
  • order::Int: The order of spatial interpolation (1 through 5 are supported).
source
FixedMeshRefinement.step!Method
step!(
    grid::Grid,
    f::Function,
    p;
    mongwane::Bool=false,
    apply_trans_zone::Bool=false,
)

Advance the solution on the entire grid by one time step of the coarsest level. This function orchestrates the time stepping across all refinement levels, handling subcycling, prolongation, restriction, and application of transition zones.

Arguments

  • grid::Grid: The FMR grid structure.
  • f::Function: The function that computes the right-hand side of the ODEs. It should have the signature f(level, k, u, p, t).
  • p: Parameters to be passed to the RHS function f.
  • mongwane::Bool: If true, enables Mongwane's subcycling method. Defaults to false.
  • apply_trans_zone::Bool: If true, applies transition zones to smooth inter-grid boundaries. Defaults to false.
source
FixedMeshRefinement.time_interpolate!Method
time_interpolate!(res, u, order)

Performs time interpolation for prolongation. u should be ordered from oldest to newest. The result is stored in res.

Arguments

  • res: The destination array for the interpolated values.
  • u: A vector of states at different time levels, ordered from oldest to newest.
  • order::Int: The order of time interpolation (1 through 4 are supported).
source
FixedMeshRefinement.transition_profileMethod
transition_profile(xl, xh, x; type=1)

Compute a transition profile value at x between xl and xh. This is used to smoothly blend solutions in transition zones.

Arguments

  • xl::Float64: The lower bound of the transition interval.
  • xh::Float64: The upper bound of the transition interval.
  • x::Float64: The point at which to evaluate the profile.
  • type::Int: The type of transition profile (1: boxstep, 2: smoothstep, 3: smootherstep).
source
FixedMeshRefinement.typetolMethod
typetol(TF::Type{<:AbstractFloat})

Return a tolerance value for a given float type TF. The tolerance is calculated as eps(TF)^(4/5). This is useful for setting default tolerances in numerical algorithms where machine precision is a factor.

Arguments

  • TF::Type{<:AbstractFloat}: The floating point type.

Returns

  • A tolerance value of the same floating point type.

Examples

julia> typetol(Float64)
3.666852862501036e-11

julia> typetol(Float32)
2.422181f-5
source