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 that refers to the mean value. Stochastic error control has to estimate which stochastic mode needs to be refined in the sense that either the polynomial degree is increased or neighbouring stochastic modes are activated. Both is represented by the multi-indices. Then unified error control allows to perform residual-based error estimation for the subresiduals that are associated to each multi-index and depend 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}}}
computes the residual-based a posteriori error estimator for the current solution vector sol
for the given model problem and stochastic coefficient C
and returns:
eta4modes
= array with total error estimators for each multi-index (corresponding to the enriched set of multi-indices with the current active modes coming first)eta4cell
= array with total error estimators for each cell in grid (for spatial refinement)multi_indices_extended
= enriched set of multi-indices used for the computations
Multi index management
Depending on the model problem and stochastic coefficient the amount of multi indices that should be added to the error estimation varies. Here are some methods that help with enriching the set of multi-indices.
ExtendableASGFEM.add_boundary_modes
— Methodadd_boundary_modes(
multi_indices;
p_extension,
tail_extension
) -> Any
adds new stochastic modes:
- for each existing mode all neighbouring modes are added (= all possible copies of that mode where one dimension is increased by 1)
p_extension
additionally ensures that the polynomial degree of the first dimension is increased by this amounttail_extension
activatestail_expansion[1]
many new stochastic modes with order 1, for each existing mode also the nexttail_expansion[2]
many higher stochastic modes are activated or increased
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 these categories:
active_int
= mode and all its direct neighbours (+1 in all active and the next dimensions) are in active_modesactive_bnd
= mode but not all its direct neighbours are in active_modesinactive_bnd
= inactive mode that has an active neighbourinactive_bnd2
= inactive mode that has an neighbour in inactive_bndinactive_else
= all other
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
There is also a hierarchical Monte carlo error estimator available that compares the solution with a higher order discrete solution for sampled deterministic problems. This is merely intended as a way to compute the 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 model problem problem
by Monte carlo sampling: for each sample a discrete finite element solution of the deterministic model problem with fixed (sampled) coefficient with polynomial order order
is computed and compared to the given stochastic Galerkin solution SolutionSGFEM
. Return values (each arrays of length M+1 where M is the length of the multi-indices) are:
totalerrorL2stress_weighted
: mean L2 error of the stress with samples weighted by the distributiontotalerrorL2u_weighted
: mean L2 error with samples weighted by the distributiontotalerrorL2stress_uniform
: mean L2 error of the stress with samples weighted uniformlytotalerrorL2u_uniform
: mean L2 error with samples weighted uniformly
The m-th component of these arrays are the errors when only multi-indices of up to order m are included.