GridapROMs.RBTransient
Reduced-basis infrastructure for transient parametric PDEs, extending RBSteady to the space–time setting.
Reduction methods
High-dimensional (transient) reduction methods extend the steady Reduction to handle the time dimension. The abstract supertype HighDimReduction dispatches to one of three strategies depending on the arguments: SteadyReduction, KroneckerReduction, and SequentialReduction.
Hyper-reduction methods
Transient hyper-reduction wraps a HighDimReduction together with a TimeCombination that encodes the ODE time-marching coefficients. See HighDimHyperReduction for details on how the θ-method (and higher-order schemes) lead to per-operator combination orders. Concrete subtypes include HighDimMDEIMHyperReduction, HighDimSOPTHyperReduction, and HighDimRBFHyperReduction.
Bases construction
The tucker function computes a Tucker decomposition of multi-dimensional snapshot arrays.
Transient projections
TransientProjection wraps spatial and temporal projections for space–time reduced-basis problems, while TransientIntegrationDomain provides the corresponding reduced integration domain.
Full API
GridapROMs.RBTransient — Module
module RBTransientReduced-basis infrastructure for transient parametric PDEs.
Extends RBSteady to the space–time setting. The additional complexity arises from the time dimension: snapshots are matrices (space × time) and bases must capture temporal as well as spatial structure. Key extensions:
Reduction methods —
TransientReduction,KroneckerReduction,SequentialReductionextend the steady reductions;HighDimReduction,SteadyReductionhandle the "high-dimensional" (full-space) and purely-spatial cases. Corresponding hyper-reduction variants:TransientHyperReduction,HighDimMDEIMHyperReduction,HighDimSOPTHyperReduction.Tucker / Kronecker bases (
BasesConstruction.jl) —tuckercomputes a Tucker decomposition, used for Kronecker-product basis representations.Transient projections (
Projections.jl) —TransientProjectionwraps a spatial projection with a temporal one.Transient integration domains (
IntegrationDomains.jl) —TransientIntegrationDomainextends DEIM-style reduced integration to include a reduced time-index set.Transient interpolations (
Interpolations.jl) —TransientGreedyInterpolation,TransientRBFInterpolation,TransientBlockInterpolation.Transient operators (
ReducedOperators.jl) —TransientRBOperatoradds time-stepping to the reduced operator interface.Time marching (
ParamTimeMarching.jl) — hooks intoParamODEsto drive the online transient solve with the reduced operator.Post-processing and I/O (
PostProcess.jl) — transient counterpart of the steady post-processing utilities.
All space-only functionality (RB spaces, hyper-reduction infrastructure, linear algebra) is imported directly from RBSteady.
GridapROMs.RBTransient.HighDimHyperReduction — Type
abstract type HighDimHyperReduction{A} <: HyperReduction{A} endHyper-reduction strategies employed in high-order (e.g. transient) problems.
Every concrete subtype stores a TimeCombination, which encodes the way an ODE time-marching scheme combines contributions from different time levels. See TimeCombination and CombinationOrder for the full treatment; here we summarise the key idea for the simplest case.
Theta method (first-order ODE)
Consider
\[M \dot{u}(t) + A\, u(t) = f(t).\]
The $\theta$-method reads
\[M \frac{u_{n+1} - u_n}{\Delta t} + \theta\, A\, u_{n+1} + (1-\theta)\, A\, u_n = f_{n+\theta},\]
which can be rewritten as
\[\left( \frac{1}{\Delta t} M + \theta\, A \right) u_{n+1} = \left( \frac{1}{\Delta t} M - (1-\theta)\, A \right) u_n + f_{n+\theta}.\]
The scheme is therefore purely one-step:
\[\boxed{ \left( \frac{1}{\Delta t} M + \theta\, A \right) u_{n+1} = \left( \frac{1}{\Delta t} M - (1-\theta)\, A \right) u_n + f_{n+\theta}. }\]
Within the ROM framework the two operators $A$ and $M$ are associated with distinct CombinationOrder indices (1 and 2 for a first-order problem). The TimeCombination object stores the scheme parameters ($\theta$, $\Delta t$, …) and the function get_coefficients returns the per-order weights that combine snapshots from successive time levels. Higher-order schemes (Newmark / Generalized-α) follow the same pattern with additional combination orders.
GridapROMs.RBTransient.HighDimMDEIMHyperReduction — Type
struct HighDimMDEIMHyperReduction{A,R<:Reduction{A,EuclideanNorm}} <: HighDimHyperReduction{A}Transient hyper-reduction based on the Matrix Discrete Empirical Interpolation Method (MDEIM). Combines a spatial HighDimReduction with a TimeCombination encoding the ODE time-marching coefficients.
Fields
reduction::R: the underlying spatial reduction.combination::TimeCombination: time-marching combination.
GridapROMs.RBTransient.HighDimRBFHyperReduction — Type
struct HighDimRBFHyperReduction{A,R<:Reduction{A,EuclideanNorm}} <: HighDimHyperReduction{A}Transient hyper-reduction based on radial basis function (RBF) interpolation. In addition to the spatial HighDimReduction and TimeCombination, it stores an AbstractRadialBasis strategy that governs the RBF kernel.
Fields
reduction::R: the underlying spatial reduction.combination::TimeCombination: time-marching combination.strategy::AbstractRadialBasis: radial basis function kernel (defaultPHS()).
GridapROMs.RBTransient.HighDimReduction — Type
abstract type HighDimReduction{A<:ReductionStyle,B<:NormStyle} <: Reduction{A,B} endAbstract supertype for reduction methods in high-order (e.g. transient) parametric problems.
Concrete subtypes:
SteadyReduction— wraps a steadyReduction; no temporal compression is applied (the ROM still time-steps).KroneckerReduction— builds a Kronecker (Tucker) product space from independent spatial and temporal reductions.SequentialReduction— uses TT-SVD (tensor-train) decomposition for the snapshot tensor.
Use the generic constructor HighDimReduction(args...; kwargs...) to dispatch to the appropriate subtype based on the arguments.
GridapROMs.RBTransient.HighDimSOPTHyperReduction — Type
struct HighDimSOPTHyperReduction{A,R<:Reduction{A,EuclideanNorm}} <: HighDimHyperReduction{A}Transient hyper-reduction based on the SOPT (Second-Order Proper Transformation) strategy. Stores a spatial HighDimReduction and a TimeCombination.
Fields
reduction::R: the underlying spatial reduction.combination::TimeCombination: time-marching combination.
GridapROMs.RBTransient.KroneckerProjection — Type
struct KroneckerProjection <: TransientProjection
projection_space::Projection
projection_time::Projection
endProjection operator for transient problems, containing a spatial projection and a temporal one. The space-time projection operator is equal to
projection_time ⊗ projection_space
which, for efficiency reasons, is never explicitly computed
GridapROMs.RBTransient.KroneckerReduction — Type
struct KroneckerReduction{A,B} <: HighDimReduction{A,B}
reductions::AbstractVector{<:Reduction}
endWrapper for reduction methods in high order problems, such as transient ones. The reduced subspaces are constructed as Kronecker product spaces
GridapROMs.RBTransient.SequentialReduction — Type
struct SequentialReduction{A,B} <: HighDimReduction{A,B}
reduction::Reduction{A,B}
endWrapper for sequential reduction methods in high-order problems, e.g. TT-SVD in transient applications
GridapROMs.RBTransient.SteadyReduction — Type
struct SteadyReduction{A,B} <: HighDimReduction{A,B}
reduction::Reduction{A,B}
endWrapper for steady reduction methods in high order problems, such as transient ones. In practice, the resulting ROM will still need to run the time marching scheme, since no temporal reduction occurs.
GridapROMs.RBTransient.TransientIntegrationDomain — Type
struct TransientIntegrationDomain{A<:TransientIntegrationDomainStyle,Ti<:Integer} <: IntegrationDomain
domain_style::A
domain_space::IntegrationDomain
indices_time::Vector{Ti}
endIntegration domain for a projection operator in a transient problem
GridapROMs.RBTransient.TransientProjection — Type
abstract type TransientProjection <: Projection endAbstract type for projection operators in transient reduced-basis problems. A TransientProjection wraps a spatial projection and a temporal projection and provides methods to project and reconstruct space–time snapshot vectors.
Interface
get_projection_space(a): return the spatialProjection.get_projection_time(a): return the temporalProjection.get_basis_space(a)/get_basis_time(a): retrieve the spatial / temporal basis matrices.
GridapROMs.RBTransient.time_enrichment — Method
time_enrichment(red::SupremizerReduction,a_primal::Projection,basis_dual;kwargs...) -> AbstractMatrixTemporal supremizer enrichment. (Approximate) Procedure:
- for every
b_dual ∈ Col(basis_dual) - compute
Φ_primal_dual = get_basis(a_primal)'*get_basis(b_dual) - compute
v = kernel(Φ_primal_dual) - compute
v′ = orth_complement(v,a_primal) - enrich
a_primal = [a_primal,v′]
GridapROMs.RBTransient.tucker — Method
tucker(red, A)
tucker(red, A, X...)Compute a Tucker decomposition of the multi-dimensional array A by sequentially applying a truncated POD (via tucker_loop) along each unfolding mode. The vector red provides a Reduction strategy for every mode (its length must equal ndims(A) - 1). When sparse matrices X are supplied they are used as inner-product weights for the corresponding modes.
Returns a Vector{Matrix{T}} of orthonormal basis matrices, one per mode.