Solver

Solving a problem requires both spatial and stochastic discretization. These are combined into a specialized vector structure, which is then passed to a solver function that executes an iterative algorithm tailored to each model problem.

SGFEVector

The spatial discretization is defined by a single finite element space from ExtendableFEM.jl, while the stochastic discretization uses a tensorized basis for the parameter space of the stochastic coefficient. Both components must be set up in advance.

Note

Currently, it is not possible to use different finite element spaces for different multi-indices. This feature may be added in the future.

ExtendableASGFEM.SGFEVectorType
struct SGFEVector{T, Tv, Ti, ONBType<:ONBasis, MIType}

A structure that extends ExtendableFEMBase.FEVector to include information about the stochastic discretization, specifically the associated TensorizedBasis.

Fields

  • FES_space: Array of finite element spaces for all stochastic modes.
  • TB: The tensorized basis used for the stochastic discretization.
  • active_modes: Indices of the active stochastic modes (with respect to TB.multi_indices).
  • length4modes: Offsets for each mode in the global vector.
  • FEVectorBlocks: Vector block mask for each multi-index (stochastic mode).
  • entries: The full coefficient vector containing all degrees of freedom.
  • last_sample: The last sample used for evaluation (for efficient repeated evaluation).
  • FEV: An FEVector used for evaluating the SGFEVector at a given sample.

This structure enables efficient storage, evaluation, and manipulation of solutions in stochastic Galerkin finite element methods.

source
ExtendableASGFEM.SGFEVectorMethod
SGFEVector(
    FES::Array{<:ExtendableFEMBase.FESpace{Tv, Ti}, 1},
    TB::TensorizedBasis{Tv, ONBType, MIType};
    active_modes,
    T,
    unames
) -> SGFEVector

Constructs an SGFEVector for the given spatial finite element spaces FES and the tensorized basis TB representing the stochastic discretization.

Arguments

  • FES: Array of finite element spaces (one for each unknown).
  • TB: The tensorized basis for the stochastic discretization.
  • active_modes: Indices of the active stochastic modes to include (default: all modes in TB).
  • T: The floating-point type for the coefficient vector (default: Tv).
  • unames: Names or identifiers for the unknowns (default: 1:length(FES)).

Returns

An SGFEVector object that stores the spatial and stochastic discretization, the coefficient vector, and all necessary metadata for efficient evaluation and manipulation in stochastic Galerkin FEM.

source
Base.getindexMethod
getindex(
    SGFEV::SGFEVector,
    u::Int64,
    i::Int64
) -> ExtendableFEMBase.FEVectorBlock

returns the FEVectorBlock for the i-th stochastic mode of the u-th unknown

source
Base.getindexMethod
getindex(
    SGFEV::SGFEVector,
    i::Int64
) -> ExtendableFEMBase.FEVectorBlock

returns the i-th stochastic mode

source
Base.lengthMethod
length(SGFEV::SGFEVector) -> Int64

returns the length of the full vector, i.e., the total number of degrees of freedom

source
Base.showMethod
show(io::IO, SGFEV::SGFEVector)

Custom show function for FEVector that prints some information on its blocks.

source
Base.sizeMethod
size(SGFEV::SGFEVector)

returns a tuple with the number of active stochastic modes and the number of spatial degrees of freedom

source
ExtendableASGFEM.set_sample!Method
set_sample!(SGFEV::SGFEVector, S::AbstractVector)

Evaluates the SGFEVector at the given sample S and stores the result in SGFEV.FEV.

Arguments

  • SGFEV: The stochastic Galerkin finite element vector to evaluate.
  • S: A vector representing the sample (values for the stochastic variables).

Details

  • The sample S is stored in SGFEV.last_sample (truncated or padded as needed).
  • The tensorized basis is evaluated at S, and the resulting coefficients are used to assemble the spatial solution in SGFEV.FEV.
  • This enables efficient evaluation of the SGFEM solution at arbitrary points in the stochastic parameter space.
source

Solve Function

ExtendableASGFEM.solve!Method
solve!(
    ::Type{ExtendableASGFEM.AbstractModelProblem},
    sol::SGFEVector,
    C::AbstractStochasticCoefficient;
    rhs,
    use_iterative_solver,
    bonus_quadorder_f,
    bonus_quadorder_a,
    kwargs...
)

Solves the specified model problem using the provided stochastic coefficient C and right-hand side rhs, writing the solution into sol.

The sol vector communicates both the spatial and stochastic discretization, as well as any initial data required for the iterative solver.

  • If use_iterative_solver (default: true) is set, an iterative solver is used. Otherwise, the full system matrix is assembled and solved directly (note: this is very slow for large systems).
  • The parameters bonus_quadorder_f (default: 0) and bonus_quadorder_a (default: 2) allow you to increase the quadrature order for terms involving the right-hand side or the stochastic coefficient, respectively.
  • Additional keyword arguments can be passed via kwargs.

If no solver is implemented for the given model problem, an error is thrown.

source