Adaptivity and error control
Residual-based a posteriori error estimation
Spatial error estimation refers to classical residual-based error estimation for the zero-th multi-index, which corresponds to the mean value. Stochastic error control determines which stochastic mode should be refined, either by increasing the polynomial degree or by activating neighboring stochastic modes. Both are represented by multi-indices. Unified error control enables residual-based error estimation for the subresiduals associated with each multi-index, depending on the model problem (see references on the main page for details).
ExtendableASGFEM.estimate
— Methodestimate(
::Type{ExtendableASGFEM.AbstractModelProblem},
sol::SGFEVector,
C::AbstractStochasticCoefficient;
kwargs...
) -> Tuple{DataType, DataType, Vector{Vector{Int64}}}
Compute the residual-based a posteriori error estimator for a stochastic Galerkin solution vector sol
.
Arguments
::Type{AbstractModelProblem}
: The model problem type for which the estimator is called. Used for dispatching to the appropriate estimator implementation.sol::SGFEVector
: The current stochastic Galerkin solution vector.C::AbstractStochasticCoefficient
: The stochastic coefficient (random field or parameterization).rhs
: (Optional) Right-hand side function for the PDE (default:nothing
).bonus_quadorder
: (Optional) Additional quadrature order for integration (default: 1).tail_extension
: (Optional) Number of additional boundary modes to include in the multi-index set (default: 5).kwargs...
: Additional keyword arguments passed to the estimator.
Returns
A tuple containing:
eta4modes::Vector{Float64}
: Total error estimator for each multi-index (stochastic mode), corresponding to the enriched set of multi-indices (with current active modes first).eta4cell::Matrix{Float64}
: Error estimator for each cell in the spatial grid (for spatial refinement), for each multi-index.multi_indices_extended
: The enriched set of multi-indices used in the computation (including boundary extensions).
Description
This function computes a residual-based a posteriori error estimator for the current SGFEM solution. It supports both spatial and stochastic adaptivity by providing error indicators for each cell and each stochastic mode. The estimator is tailored to the model problem and the stochastic coefficient, and can be extended to include additional boundary modes for improved reliability.
If no specialized estimator is available for the given model problem type, an error is raised and empty arrays are returned. ```
Multi-index management
Depending on the model problem and stochastic coefficient, the number of multi-indices to be added for error estimation varies. The following methods assist in enriching the set of multi-indices.
ExtendableASGFEM.add_boundary_modes
— Methodadd_boundary_modes(
multi_indices;
p_extension,
tail_extension
) -> Any
Adds new stochastic modes (multi-indices) to the current set by extending in several ways:
- For each existing mode, all neighboring modes are added (i.e., all possible copies of that mode where one dimension is increased by 1).
p_extension
ensures that the polynomial degree in the first dimension is increased by this amount (adds new modes with higher degree in the first component).tail_extension
is a two-element vector:tail_extension[1]
: Adds this many new modes of order 1 (i.e., with a single nonzero entry).tail_extension[2]
: For each existing mode, also activates or increases the nexttail_extension[2]
higher stochastic modes in each dimension.
Returns the extended set of multi-indices.
Arguments
multi_indices
: The current set of multi-indices (vector of integer vectors).p_extension
: Number of additional degrees to add in the first dimension (default: 1).tail_extension
: Two-element vector controlling the number and type of new modes to add (default:[10, 2]
).
Returns
A new vector containing the extended set of multi-indices.
ExtendableASGFEM.classify_modes
— Functionclassify_modes(multi_indices) -> NTuple{5, Vector{Int64}}
classify_modes(
multi_indices,
active_modes
) -> NTuple{5, Vector{Int64}}
Classifies all multi-indices into the following categories based on their neighborhood relations:
active_int
: Active mode where all direct neighbors (obtained by increasing any active or next dimension by 1) are also inactive_modes
.active_bnd
: Active mode where not all direct neighbors are inactive_modes
(i.e., it lies on the active boundary).inactive_bnd
: Inactive mode that has at least one active neighbor (first layer of inactive modes).inactive_bnd2
: Inactive mode that has a neighbor ininactive_bnd
(second layer of inactive modes).inactive_else
: All other inactive modes.
Arguments
multi_indices
: The set of all multi-indices to classify (vector of integer vectors).active_modes
: The set of currently active modes (default: allmulti_indices
).
Returns
A tuple of five vectors containing the indices of the multi-indices in the following order: (inactive_else, inactive_bnd, inactive_bnd2, active_bnd, active_int)
.
ExtendableASGFEM.generate_multiindices
— Methodgenerate_multiindices(M, deg) -> Any
generates multi-indices for M
dimensions up to degree deg
ExtendableASGFEM.prepare_multi_indices!
— Methodprepare_multi_indices!(multi_indices; minimal_length)
ensures that all multi-indices have the same length
Monte Carlo sampling estimator
A hierarchical Monte Carlo error estimator is also available. It compares the solution with a higher-order discrete solution for sampled deterministic problems. This is mainly intended to compute a reference error to assess the efficiency of the residual-based error estimator.
ExtendableASGFEM.calculate_sampling_error
— Methodcalculate_sampling_error(
SolutionSGFEM::SGFEVector,
C::AbstractStochasticCoefficient;
problem,
bonus_quadorder_a,
bonus_quadorder_f,
order,
rhs,
dim,
Msamples,
parallel_sampling,
dimensionwise_error,
energy_norm,
debug,
nsamples
) -> NTuple{4, Vector{Float64}}
Estimates the error for the given model problem by Monte Carlo sampling. For each sample, a deterministic finite element solution is computed for a fixed (sampled) coefficient with polynomial order order
and compared to the provided stochastic Galerkin solution SolutionSGFEM
.
Returns
A tuple of arrays (each of length M+1
, where M
is the maximum order of the multi-indices):
totalerrorL2stress_weighted
: Mean L2 error of the stress, weighted by the distribution.totalerrorL2u_weighted
: Mean L2 error of the solution, weighted by the distribution.totalerrorL2stress_uniform
: Mean L2 error of the stress, weighted uniformly.totalerrorL2u_uniform
: Mean L2 error of the solution, weighted uniformly.
The m-th component of these arrays gives the error when only multi-indices up to order m are included.
Arguments
SolutionSGFEM
: The stochastic Galerkin solution (SGFEVector).C
: The stochastic coefficient.problem
: The model problem type (default:LogTransformedPoissonProblemPrimal
).bonus_quadorder_a
: Additional quadrature order for the coefficient (default: 10).bonus_quadorder_f
: Additional quadrature order for the right-hand side (default: 0).order
: Polynomial order for the deterministic reference solution (default: 2).rhs
: Right-hand side function (optional).dim
: Spatial dimension (default: inferred from the solution).Msamples
: Number of random variables to sample (default:maxm(C)
).parallel_sampling
: Whether to use parallel sampling (default: true).dimensionwise_error
: If true, computes errors dimensionwise (default: false).energy_norm
: If true, uses the energy norm for stress error (default: true).debug
: If true, enables debug output and plotting (default: false).nsamples
: Number of Monte Carlo samples (default: 100).
Details
This function computes deterministic reference solutions for each sample, then compares them to the SGFEM solution to estimate the error. Both weighted and unweighted (uniform) averages are returned for stress and solution errors.