FixedMeshRefinement
Documentation for FixedMeshRefinement.
FixedMeshRefinement.Grid
FixedMeshRefinement.Level
Base.show
Base.show
Base.show
Base.show
FixedMeshRefinement.apply_transition_zone!
FixedMeshRefinement.calc_Yn_from_kcs!
FixedMeshRefinement.cycle_state!
FixedMeshRefinement.durationstring
FixedMeshRefinement.fill_buffer!
FixedMeshRefinement.fine_to_coarse_index
FixedMeshRefinement.get_boundary_indices
FixedMeshRefinement.get_diagnostic_state
FixedMeshRefinement.get_interior_indices
FixedMeshRefinement.get_maximum_grid_points
FixedMeshRefinement.get_offset_indices
FixedMeshRefinement.get_rhs_evaluation_indices
FixedMeshRefinement.get_rk4_evaluation_indices
FixedMeshRefinement.get_rk4_tmp_state
FixedMeshRefinement.get_rk_stage
FixedMeshRefinement.get_state
FixedMeshRefinement.get_tmp_state
FixedMeshRefinement.get_total_grid_points
FixedMeshRefinement.get_x
FixedMeshRefinement.integrate_simpson
FixedMeshRefinement.integrate_trapezoidal
FixedMeshRefinement.merge_grid_levels
FixedMeshRefinement.prolongate!
FixedMeshRefinement.prolongate_mongwane!
FixedMeshRefinement.restrict_injection!
FixedMeshRefinement.rk4!
FixedMeshRefinement.rk4_dense_output_y!
FixedMeshRefinement.searchsortednearest
FixedMeshRefinement.shift_grid_boundaries!
FixedMeshRefinement.shift_level_boundaries!
FixedMeshRefinement.spatial_interpolate!
FixedMeshRefinement.step!
FixedMeshRefinement.time_interpolate!
FixedMeshRefinement.transition_profile
FixedMeshRefinement.typetol
FixedMeshRefinement.Grid
— TypeGrid
A struct representing the entire FMR grid, which consists of multiple Level
s.
Fields
levels::Vector{Level}
: A vector ofLevel
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.
FixedMeshRefinement.Level
— TypeLevel
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 forOffsetArray
s of this level.
Base.show
— MethodBase.show(io::IO, grid::Grid)
Display a compact summary of the Grid
.
Base.show
— MethodBase.show(io::IO, level::Level)
Display a compact summary of the Level
.
Base.show
— MethodBase.show(io::IO, ::MIME"text/plain", grid::Grid)
Display a detailed summary of the Grid
structure.
Base.show
— MethodBase.show(io::IO, ::MIME"text/plain", level::Level)
Display a detailed summary of the Level
.
FixedMeshRefinement.apply_transition_zone!
— Methodapply_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.
FixedMeshRefinement.calc_Yn_from_kcs!
— Methodcalc_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 calculatedY
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.
FixedMeshRefinement.cycle_state!
— Methodcycle_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.
FixedMeshRefinement.durationstring
— Methoddurationstring(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
FixedMeshRefinement.fill_buffer!
— Methodfill_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
: Iftrue
, fill the buffer withNaN
to avoid accidental reuse. Defaults totrue
.
FixedMeshRefinement.fine_to_coarse_index
— Methodfine_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.
FixedMeshRefinement.get_boundary_indices
— Methodget_boundary_indices(level::Level) -> NTuple{2,StepRange{Int,Int}}
Return the indices of the boundary points on each side (ghost or buffer).
FixedMeshRefinement.get_diagnostic_state
— Methodget_diagnostic_state(level::Level)
Return the diagnostic state variables of the level
as an OffsetArray
.
FixedMeshRefinement.get_interior_indices
— Methodget_interior_indices(level::Level) -> UnitRange{Int}
Return the indices of the interior grid points.
FixedMeshRefinement.get_maximum_grid_points
— Methodget_maximum_grid_points(level::Level) -> Int
Return the maximum number of grid points at this level (including ghost, buffer and unused points).
FixedMeshRefinement.get_offset_indices
— Methodget_offset_indices(level::Level) -> UnitRange{Int}
Return the indices of the grid points in the OffsetArray
of the level
.
FixedMeshRefinement.get_rhs_evaluation_indices
— Methodget_rhs_evaluation_indices(level::Level) -> UnitRange{Int}
Return the indices that requires evaluation of the right-hand side.
FixedMeshRefinement.get_rk4_evaluation_indices
— Methodget_rk4_evaluation_indices(level::Level) -> UnitRange{Int}
Return the indices that requires evaluation of the RK4 stages.
FixedMeshRefinement.get_rk4_tmp_state
— Methodget_rk4_tmp_state(level::Level)
Return a temporary array with the same size as the state variables of the level
as an OffsetArray
.
FixedMeshRefinement.get_rk_stage
— Methodget_rk_stage(level::Level, i::Int)
Return the i
-th intermediate state k_i
for the Runge-Kutta time stepping scheme as an OffsetArray
.
FixedMeshRefinement.get_state
— Functionget_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.
FixedMeshRefinement.get_tmp_state
— Methodget_tmp_state(level::Level)
Return the temporary state variables of the level
as an OffsetArray
.
FixedMeshRefinement.get_total_grid_points
— Methodget_total_grid_points(level::Level) -> Int
Return the total number of grid points at this level (including ghost and buffer points).
FixedMeshRefinement.get_x
— Methodget_x(level::Level)
Return the grid point coordinates of the level
as an OffsetArray
.
FixedMeshRefinement.integrate_simpson
— Methodintegrate_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
.
FixedMeshRefinement.integrate_trapezoidal
— Methodintegrate_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
.
FixedMeshRefinement.merge_grid_levels
— Methodmerge_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.
FixedMeshRefinement.prolongate!
— Methodprolongate!(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.
FixedMeshRefinement.prolongate_mongwane!
— Methodprolongate_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.
FixedMeshRefinement.restrict_injection!
— Methodrestrict_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
: Iftrue
, the restriction is only applied outside of the transition zones. Defaults tofalse
.
FixedMeshRefinement.rk4!
— Methodrk4!(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 signaturef(level, k, u, p, t)
.p
: Parameters to be passed to the RHS functionf
.mongwane::Bool
: Iftrue
, enables special buffer filling for Mongwane's subcycling method. Defaults tofalse
.
FixedMeshRefinement.rk4_dense_output_y!
— Methodrk4_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 stagesk1
throughk4
.
FixedMeshRefinement.searchsortednearest
— Methodsearchsortednearest(a::AbstractVector, x)
Find the index of the element in a sorted array a
that is nearest to x
.
FixedMeshRefinement.shift_grid_boundaries!
— Methodshift_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 ofnum_shift_points::NTuple{2,Int}
: Number of points to shift on left and right boundariesfunc_fill_extended::Function=fill_extended_grid_extrapolate!
: Function to compute fill values, takes index and level as arguments
FixedMeshRefinement.shift_level_boundaries!
— Methodshift_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.
FixedMeshRefinement.spatial_interpolate!
— Methodspatial_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 inu
to center the interpolation stencil.order::Int
: The order of spatial interpolation (1 through 5 are supported).
FixedMeshRefinement.step!
— Methodstep!(
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 signaturef(level, k, u, p, t)
.p
: Parameters to be passed to the RHS functionf
.mongwane::Bool
: Iftrue
, enables Mongwane's subcycling method. Defaults tofalse
.apply_trans_zone::Bool
: Iftrue
, applies transition zones to smooth inter-grid boundaries. Defaults tofalse
.
FixedMeshRefinement.time_interpolate!
— Methodtime_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).
FixedMeshRefinement.transition_profile
— Methodtransition_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).
FixedMeshRefinement.typetol
— Methodtypetol(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