GridapROMs.RBSteady

GridapROMs.RBSteady.AbstractMDEIMReductionType
abstract type AbstractMDEIMReduction{A} <: Reduction{A,EuclideanNorm} end

Type representing a hyper-reduction approximation by means of a MDEIM algorithm. Check this for more details on MDEIM. This reduction strategy is usually employed only for the approximation of a residual and/or Jacobian of a differential problem. Note that orthogonality with respect to a norm other than the euclidean is not required for this reduction type.

Subtypes:

source
GridapROMs.RBSteady.AdaptiveReductionType
struct AdaptiveReduction{A,B,R<:DirectReduction{A,B}} <: GreedyReduction{A,B}
  reduction::R
  adaptive_nparams::Int
  adaptive_tol::Float64
  adaptive_maxiter::Int
end

Not implemented yet. Will serve as a parameter-adaptivity greedy reduction algorithm

source
GridapROMs.RBSteady.AffineReductionType
struct AffineReduction{A,B} <: DirectReduction{A,B}
  red_style::A
  norm_style::B
end

Reduction employed when the input data is independent with respect to the considered realization. Therefore, simply considering a number of parameters equal to 1 suffices for this type of reduction

source
GridapROMs.RBSteady.DofMapCoreType
struct DofMapCore{T,A<:AbstractArray{T,3},B<:AbstractArray} <: AbstractTTCore{T,3}
  core::A
  dof_map::B
end

Represents a tensor train core core reindexed by means of an index mapping dof_map

source
GridapROMs.RBSteady.EnergyNormType
struct EnergyNorm <: NormStyle
  norm_op::Function
end

Trait indicating that the reduction algorithm will produce a basis orthogonal in the norm specified by norm_op. Note: norm_op should represent a symmetric, positive definite bilinear form (matrix)

source
GridapROMs.RBSteady.EvalRBSpaceType
struct EvalRBSpace{A<:RBSpace,B<:AbstractRealization} <: RBSpace
  subspace::A
  realization::B
end

Conceptually this isn't needed, but it helps dispatching according to steady/transient cases

source
GridapROMs.RBSteady.FixedSVDRankType
struct FixedSVDRank <: ReductionStyle
  rank::Int
end

Struct employed when the chosen reduction algorithm is a truncated POD at a rank rank. Check this reference for more details on the truncated POD algorithm

source
GridapROMs.RBSteady.GenericRBOperatorType
struct GenericRBOperator{O,A} <: RBOperator{O}
  op::ParamOperator{O}
  trial::RBSpace
  test::RBSpace
  lhs::A
  rhs::AffineContribution
end

Fields:

  • op: underlying high dimensional FE operator
  • trial: reduced trial space
  • test: reduced trial space
  • lhs: hyper-reduced left hand side
  • rhs: hyper-reduced right hand side
source
GridapROMs.RBSteady.HRParamArrayType
struct HRParamArray{T,N,A,B,C<:ParamArray{T,N}} <: ParamArray{T,N}
  fecache::A
  coeff::B
  hypred::C
end

Parametric vector returned after the online phase of a hyper-reduction strategy. Fields:

  • fecache: represents a parametric residual/Jacobian computed via integration on an IntegrationDomain

  • coeff: parameter-dependent coefficient computed during the online phase according to the formula

    coeff = Φi⁻¹ fecache[i,:]

    where (Φi,i) are the interpolation and the reduced integration domain stored in a HyperReduction object.

  • hypred: the ouptut of the online phase of a hyper-reduction strategy, acoording to the formula

    hypred = Φrb * coeff

    where Φrb is the basis stored in a HyperReduction object

source
GridapROMs.RBSteady.HyperReductionType
abstract type HyperReduction{
  A<:Reduction,
  B<:ReducedProjection,
  C<:IntegrationDomain
  } <: Projection end

Subtype of a Projection dedicated to the outputd of a hyper-reduction (e.g. an empirical interpolation method (EIM)) procedure applied on residual/jacobians of a differential problem. This procedure can be summarized in the following steps:

  1. compute a snapshots tensor T
  2. construct a Projection Φ by running the function reduction on T
  3. find the EIM quantities (Φi,i), by running the function empirical_interpolation

on Φ

The triplet (Φ,Φi,i) represents the minimum information needed to run the online phase of the hyper-reduction. However, we recall that a RB method requires the (Petrov-)Galerkin projection of residuals/Jacobianson a reduced subspace built from solution snapshots, instead of providing the projection Φ we return the reduced projection Φrb, where

  • for residuals: Φrb = test_basisᵀ Φ
  • for Jacobians: Φrb = test_basisᵀ Φ trial_basis

The output of this operation is a ReducedProjection. Therefore, a HyperReduction is completely characterized by the triplet (Φrb,Φi,i). Subtypes:

source
GridapROMs.RBSteady.LRApproxRankType
struct LRApproxRank <: ReductionStyle
  opts::LRAOptions
end

Struct employed when the chosen reduction algorithm is a randomized POD that leverages the package LowRankApprox. The field opts specifies the options needed to run the randomized POD

source
GridapROMs.RBSteady.MatrixDomainType
struct MatrixDomain{T,S} <: IntegrationDomain
  cells::Vector{Int32}
  cell_irows::Table{T,Vector{T},Vector{Int32}}
  cell_icols::Table{S,Vector{S},Vector{Int32}}
end

Integration domain for a projection vector operator in a steady problem

source
GridapROMs.RBSteady.NormedProjectionType
struct NormedProjection <: Projection
  projection::Projection
  norm_matrix::MatrixOrTensor
end

Represents a Projection projection spanning a space equipped with an inner product represented by the quantity norm_matrix

source
GridapROMs.RBSteady.PODReductionType
struct PODReduction{A,B} <: DirectReduction{A,B}
  red_style::A
  norm_style::B
  nparams::Int
end

Reduction by means of a truncated POD. The field nparams indicates the number of parameters selected for the computation of the snapshots

source
GridapROMs.RBSteady.RBParamVectorType
struct RBParamVector{T,A<:ParamVector{T},B} <: ParamArray{T,1}
  data::A
  fe_data::B
end

Parametric vector obtained by applying a Projection on a high-dimensional parametric FE vector fe_data, which is stored (but mostly unused) for conveniency

source
GridapROMs.RBSteady.RBSolverType
struct RBSolver{A<:GridapType,B} <: GridapType
  fesolver::A
  state_reduction::Reduction
  residual_reduction::Reduction
  jacobian_reduction::B
end

Wrapper around a FE solver (e.g. NonlinearSolver or ODESolver in Gridap) with additional information on the reduced basis (RB) method employed to solve a given problem dependent on a set of parameters. A RB method is a projection-based reduced order model where

  1. a suitable subspace of a FESpace is sought, of dimension n << Nₕ
  2. a matrix-based discrete empirical interpolation method (MDEIM) is performed

to approximate the manifold of the parametric residuals and jacobians

  1. the EIM approximations are compressed with (Petrov-)Galerkin projections

onto the subspace

  1. for every desired choice of parameters, numerical integration is performed, and

the resulting n × n system of equations is cheaply solved

In particular:

  • tol: tolerance used in the projection-based truncated proper orthogonal decomposition (TPOD) or in the tensor train singular value decomposition (TT-SVD), where a basis spanning the reduced subspace is computed; the value of tol is responsible for selecting the dimension of the subspace, i.e. n = n(tol)
  • nparams_state: number of snapshots considered when running TPOD or TT-SVD
  • nparams_res: number of snapshots considered when running MDEIM for the residual
  • nparams_jac: number of snapshots considered when running MDEIM for the jacobian
  • nparams_test: number of snapshots considered when computing the error the RB method commits with respect to the FE procedure
source
GridapROMs.RBSteady.RBVectorType
struct RBVector{T,A<:AbstractVector{T},B} <: AbstractVector{T}
  data::A
  fe_data::B
end

Vector obtained by applying a Projection on a high-dimensional FE vector fe_data, which is stored (but mostly unused) for conveniency

source
GridapROMs.RBSteady.ROMPerformanceType
struct ROMPerformance
  error
  speedup
end

Allows to compute errors and computational speedups to compare the properties of the algorithm with the FE performance.

source
GridapROMs.RBSteady.SearchSVDRankType
struct SearchSVDRank <: ReductionStyle
  tol::Float64
end

Struct employed when the chosen reduction algorithm is a truncated POD at a tolerance tol. Check this reference for more details on the truncated POD algorithm

source
GridapROMs.RBSteady.SparseCoreCSCType
struct SparseCoreCSC{T,Ti} <: SparseCore{T,3}
  array::Array{T,3}
  sparsity::SparsityCSC{T,Ti}
end

Tensor train cores for sparse matrices in CSC format

source
GridapROMs.RBSteady.SparseCoreCSC4DType
struct SparseCoreCSC4D{T,Ti} <: SparseCore{T,4}
  core::SparseCoreCSC{T,Ti}
  sparse_indexes::Vector{CartesianIndex{2}}
end

Tensor train cores for sparse matrices in CSC format, reshaped as 4D arrays

source
GridapROMs.RBSteady.SupremizerReductionType
struct SupremizerReduction{A,R<:Reduction{A,EnergyNorm}} <: Reduction{A,EnergyNorm}
  reduction::R
  supr_op::Function
  supr_tol::Float64
end

Wrapper for reduction methods reduction that require an additional step of stabilization, by means of a supremizer enrichment. Check this for more details in a steady setting, and this for more details in a transient setting. The fields supr_op and supr_tol (which is only needed only in transient applications) are respectively the supremizing operator and the tolerance involved in the enrichment. For a saddle point problem with a Jacobian of the form

[ A Bᵀ B 0 ]

this operator is given by the bilinear form representing the matrix Bᵀ.

source
GridapROMs.RBSteady.TTSVDProjectionType
struct TTSVDProjection <: Projection
  cores::AbstractVector{<:AbstractArray{T,3} where T}
  dof_map::AbstractDofMap
end

Projection stemming from a tensor train SVD ttsvd. For reindexing purposes a field dof_map is provided along with the tensor train cores cores

source
GridapROMs.RBSteady.TTSVDRanksType
struct TTSVDRanks{T<:TTSVDStyle} <: ReductionStyle
  style::Vector{<:ReductionStyle}
  unsafe::T
end

Struct employed when the chosen reduction algorithm is a TTSVD, with reduction algorithm at each step specified in the vector of reduction styles style. Check this reference for more details on the TTSVD algorithm

source
GridapROMs.RBSteady.TTSVDReductionType
struct TTSVDReduction{B} <: DirectReduction{TTSVDRanks,B}
  red_style::TTSVDRanks
  norm_style::B
  nparams::Int
end

Reduction by means of a TTSVD. The field nparams indicates the number of parameters selected for the computation of the snapshots

source
GridapROMs.RBSteady.VectorDomainType
struct VectorDomain{T} <: IntegrationDomain
  cells::Vector{Int32}
  cell_irows::Table{T,Vector{T},Vector{Int32}}
end

Integration domain for a projection vector operator in a steady problem

source
GridapROMs.RBSteady.contractionMethod
contraction(Φₗₖ::AbstractArray{T,3},Aₖ::AbstractArray{T,3}) -> AbstractArray{T,4}
contraction(Φₗₖ::AbstractArray{T,3},Aₖ::AbstractArray{T,3},Φᵣₖ::AbstractArray{T,3}) -> AbstractArray{T,6}

Contraction of tensor train cores, as a result of a (Petrov-)Galerkin projection. The contraction of Aₖ by Φₗₖ is the left-contraction of a TT core Aₖ by a (left, test) TT core Φₗₖ, whereas the contraction of Aₖ by Φᵣₖ is the right-contraction of a TT core Aₖ by a (right, trial) TT core Φᵣₖ. The dimension of the output of a contraction involving N factors is: 3N - N = 2N.

source
GridapROMs.RBSteady.contractionMethod
contraction(basis::AbstractArray,coefficient::AbstractArray) -> AbstractArray

Multiplies a reduced basis basis by a set of reduced coeffiecients coefficient. It acts as a generalized linear combination, since basis might have a dimension higher than 2.

source
GridapROMs.RBSteady.cores2basisMethod
cores2basis(cores::AbstractArray{T,3}...) -> AbstractMatrix
cores2basis(dof_map::AbstractDofMap{D},cores::AbstractArray{T,3}...) -> AbstractMatrix

Returns a basis in a matrix format from a list of tensor train cores cores. When also supplying the indexing strategy dof_map, the result is reindexed accordingly

source
GridapROMs.RBSteady.empirical_interpolationMethod
empirical_interpolation(a::Projection) -> (AbstractVector,AbstractMatrix)

Computes the EIM of a. The outputs are:

  • a vector of integers i, corresponding to a list of interpolation row indices
  • a matrix Φi = view(Φ,i), where Φ = get_basis(a). This quantity represents the restricted basis on the set of interpolation rows i
source
GridapROMs.RBSteady.enrich!Method
enrich!(
  red::SupremizerReduction,
  a::BlockProjection,
  norm_matrix::MatrixOrTensor,
  supr_matrix::MatrixOrTensor) -> Nothing

In-place augmentation of the primal block of a BlockProjection a. This function has the purpose of stabilizing the reduced equations stemming from a saddle point problem

source
GridapROMs.RBSteady.eval_performanceMethod
eval_performance(
  solver::RBSolver,
  feop::ParamOperator,
  fesnaps::AbstractSnapshots,
  rbsnaps::AbstractSnapshots,
  festats::CostTracker,
  rbstats::CostTracker
  ) -> ROMPerformance

Arguments:

  • solver: solver for the reduced problem
  • feop: FE operator representing the PDE
  • fesnaps: online snapshots of the FE solution
  • rbsnaps: reduced approximation of fesnaps
  • festats: time and memory consumption needed to compute fesnaps
  • rbstats: time and memory consumption needed to compute rbsnaps

Returns the performance of the reduced algorithm, in terms of the (relative) error between rbsnaps and fesnaps, and the computational speedup between rbstats and festats

source
GridapROMs.RBSteady.galerkin_projectionMethod
galerkin_projection(Φₗ,A) -> Any
galerkin_projection(Φₗ,A,Φᵣ) -> Any

Galerkin projection of A on the subspaces specified by a (left, test) subspace Φₗ (row projection) and a (right, trial) subspace Φᵣ (column projection)

source
GridapROMs.RBSteady.galerkin_projectionMethod
galerkin_projection(a::Projection,b::Projection) -> ReducedProjection
galerkin_projection(a::Projection,b::Projection,c::Projection,args...) -> ReducedProjection

(Petrov) Galerkin projection of a projection map b onto the subspace a (row projection) and, if applicable, onto the subspace c (column projection)

source
GridapROMs.RBSteady.gram_schmidtMethod
gram_schmidt(A::AbstractMatrix,args...) -> AbstractMatrix
gram_schmidt(A::AbstractMatrix,X::AbstractSparseMatrix,args...) -> AbstractMatrix

Gram-Schmidt orthogonalization for a matrix A under a Euclidean norm. A (positive definite) sparse matrix X representing an inner product on the row space of A can be provided to make the result orthogonal under a different norm

source
GridapROMs.RBSteady.inv_project!Method
inv_project!(x::AbstractArray,a::Projection,x̂::AbstractArray) -> Nothing

In-place recasting of a low-dimensional object the high-dimensional space in which a is immersed

source
GridapROMs.RBSteady.inv_projectMethod
inv_project(a::Projection,x::AbstractArray) -> AbstractArray

Recasts a low-dimensional object x onto the high-dimensional space in which a is immersed

source
GridapROMs.RBSteady.jacobian_snapshotsMethod
jacobian_snapshots(solver::RBSolver,op::ParamOperator,s::AbstractSnapshots;nparams) -> Contribution
jacobian_snapshots(solver::RBSolver,op::ODEParamOperator,s::AbstractSnapshots;nparams) -> Tuple{Vararg{Contribution}}

Returns a Jacobian Contribution relative to the FE operator op. The quantity s denotes the solution snapshots in which we evaluate the jacobian. Note that we can select a smaller number of parameters nparams compared to the number of parameters used to compute s. In transient settings, the output is a tuple whose nth element is the Jacobian relative to the nth temporal derivative

source
GridapROMs.RBSteady.load_operatorMethod
load_operator(dir,feop::ParamOperator;kwargs...) -> RBOperator

Given a FE operator feop, load its reduced counterpart stored in the directory dir. Throws an error if the reduced operator has not been previously saved to file

source
GridapROMs.RBSteady.load_snapshotsMethod
load_snapshots(dir;label="") -> AbstractSnapshots

Load the snapshots at the directory dir. Throws an error if the snapshots have not been previously saved to file

source
GridapROMs.RBSteady.orth_complement!Method
orth_complement!(v::AbstractVector,basis::AbstractMatrix,args...) -> Nothing

In-place orthogonal complement of v on the column space of basis. When a symmetric, positive definite matrix X is provided as an argument, the output is X-orthogonal, otherwise it is ℓ²-orthogonal

source
GridapROMs.RBSteady.orth_projectionMethod
orth_projection(v::AbstractVector, basis::AbstractMatrix, args...) -> AbstractVector

Orthogonal projection of v on the column space of basis. When a symmetric, positive definite matrix X is provided as an argument, the output is X-orthogonal, otherwise it is ℓ²-orthogonal

source
GridapROMs.RBSteady.plot_a_solutionMethod
plot_a_solution(
  dir::String,
  feop::ParamOperator,
  sol::AbstractSnapshots,
  sol_approx::AbstractSnapshots,
  args...;
  kwargs...
  ) -> Nothing

Plots a single FE solution, RB solution, and the point-wise error between the two, by selecting the first FE snapshot in sol and the first reduced snapshot in sol_approx

source
GridapROMs.RBSteady.project!Method
project!(x̂::AbstractArray,a::Projection,x::AbstractArray,args...) -> Nothing

In-place projection of a high-dimensional object x onto the subspace represented by a

source
GridapROMs.RBSteady.projectMethod
project(a::Projection,x::AbstractArray,args...) -> AbstractArray

Projects a high-dimensional object x onto the subspace represented by a

source
GridapROMs.RBSteady.projectionMethod
projection(red::Reduction,s::AbstractArray) -> Projection
projection(red::Reduction,s::AbstractArray,X::MatrixOrTensor) -> Projection

Constructs a Projection from a collection of snapshots s. An inner product represented by the quantity X can be provided, in which case the resulting Projection will be X-orthogonal

source
GridapROMs.RBSteady.reduced_jacobianMethod
reduced_jacobian(
  solver::RBSolver,
  op::ParamOperator,
  red_trial::RBSpace,
  red_test::RBSpace,
  s::AbstractSnapshots
  ) -> Union{AffineContribution,TupOfAffineContribution}

Reduces the Jacobian contained in op via hyper-reduction. This function first builds the Jacobian snapshots, which are then reduced according to the strategy reduced_jacobian specified in the reduced solver solver. In transient applications, the output is a tuple of length equal to the number of Jacobians(i.e., equal to the order of the ODE plus one)

source
GridapROMs.RBSteady.reduced_operatorMethod
reduced_operator(solver::RBSolver,feop::ParamOperator,args...;kwargs...) -> RBOperator
reduced_operator(solver::RBSolver,feop::TransientParamOperator,args...;kwargs...) -> TransientRBOperator

Computes a RB operator from the FE operator feop

source
GridapROMs.RBSteady.reduced_residualMethod
reduced_residual(
  solver::RBSolver,
  op::ParamOperator,
  red_test::RBSpace,
  s::AbstractSnapshots
  ) -> AffineContribution

Reduces the residual contained in op via hyper-reduction. This function first builds the residual snapshots, which are then reduced according to the strategy residual_reduction specified in the reduced solver solver

source
GridapROMs.RBSteady.reduced_spacesMethod
reduced_spaces(solver::RBSolver,feop::ParamOperator,s::AbstractSnapshots
  ) -> (RBSpace, RBSpace)

Computes the subspace of the test, trial FESpaces contained in the FE operator feop by compressing the snapshots s

source
GridapROMs.RBSteady.reduced_weak_formMethod
reduced_weak_form(
  solver::RBSolver,
  op::ParamOperator,
  red_trial::RBSpace,
  red_test::RBSpace,
  s::AbstractSnapshots
  ) -> (AffineContribution,Union{AffineContribution,TupOfAffineContribution})

Reduces the residual/Jacobian contained in op via hyper-reduction. Check the functions reduced_residual and reduced_jacobian for more details

source
GridapROMs.RBSteady.reductionMethod
reduction(red::Reduction,A::AbstractArray,args...) -> AbstractArray
reduction(red::Reduction,A::AbstractArray,X::AbstractSparseMatrix) -> AbstractArray

Given an array (of snapshots) A, returns a reduced basis obtained by means of the reduction strategy red

source
GridapROMs.RBSteady.residual_snapshotsMethod
residual_snapshots(solver::RBSolver,op::ParamOperator,s::AbstractSnapshots;nparams) -> Contribution
residual_snapshots(solver::RBSolver,op::ODEParamOperator,s::AbstractSnapshots;nparams) -> Contribution

Returns a residual Contribution relative to the FE operator op. The quantity s denotes the solution snapshots in which we evaluate the residual. Note that we can select a smaller number of parameters nparams compared to the number of parameters used to compute s

source
GridapROMs.RBSteady.sequential_productMethod
sequential_product(a::AbstractArray,b::AbstractArray...) -> AbstractArray

This function sequentially multiplies the results of several (sequential as well) calls to contraction

source
GridapROMs.RBSteady.solution_snapshotsMethod
solution_snapshots(solver::NonlinearSolver,feop::ParamOperator,r::Realization) -> SteadySnapshots
solution_snapshots(solver::ODESolver,feop::TransientParamOperator,r::TransientRealization,u0) -> TransientSnapshots

The problem encoded in the FE operator feop is solved several times, and the solution snapshots are returned along with the information related to the computational cost of the FE method. In transient settings, an initial condition u0 should be provided.

source
GridapROMs.RBSteady.tpodMethod
tpod(red_style::ReductionStyle,A::AbstractMatrix) -> AbstractMatrix
tpod(red_style::ReductionStyle,A::AbstractMatrix,X::AbstractSparseMatrix) -> AbstractMatrix

Truncated proper orthogonal decomposition of A. When provided, X is a (symmetric, positive definite) norm matrix with respect to which the output is made orthogonal. If X is not provided, the output is orthogonal with respect to the euclidean norm

source
GridapROMs.RBSteady.ttsvdMethod
ttsvd(red_style::TTSVDRanks,A::AbstractArray) -> AbstractVector{<:AbstractArray}
ttsvd(red_style::TTSVDRanks,A::AbstractArray,X::AbstractRankTensor) -> AbstractVector{<:AbstractArray}

Tensor train SVD of A. When provided, X is a norm tensor (representing a symmetric, positive definite matrix) with respect to which the output is made orthogonal. Note: if ndims(A) = N, the length of the ouptput is N-1, since we are not interested in reducing the axis of the parameters. Check this reference for more details

source
GridapROMs.RBSteady.union_basesMethod
union_bases(a::Projection,b::Projection,args...) -> Projection

Computes the projection corresponding to the union of a and b. In essence this operation performs as

gram_schmidt(union(get_basis(a),get_basis(b)))

source