GridapROMs.RBTransient

GridapROMs.RBTransientModule
module RBTransient

Reduced-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 methodsTransientReduction, KroneckerReduction, SequentialReduction extend the steady reductions; HighDimReduction, SteadyReduction handle the "high-dimensional" (full-space) and purely-spatial cases. Corresponding hyper-reduction variants: TransientHyperReduction, HighDimMDEIMHyperReduction, HighDimSOPTHyperReduction.

  • Tucker / Kronecker bases (BasesConstruction.jl) — tucker computes a Tucker decomposition, used for Kronecker-product basis representations.

  • Transient projections (Projections.jl) — TransientProjection wraps a spatial projection with a temporal one.

  • Transient integration domains (IntegrationDomains.jl) — TransientIntegrationDomain extends DEIM-style reduced integration to include a reduced time-index set.

  • Transient interpolations (Interpolations.jl) — TransientGreedyInterpolation, TransientRBFInterpolation, TransientBlockInterpolation.

  • Transient operators (ReducedOperators.jl) — TransientRBOperator adds time-stepping to the reduced operator interface.

  • Time marching (ParamTimeMarching.jl) — hooks into ParamODEs to 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.

source
GridapROMs.RBTransient.HighDimHyperReductionType
abstract type HighDimHyperReduction{A} <: HyperReduction{A} end

Hyper reduction strategies employed in high-order (e.g. transient) problems. They feature a field combine, a TimeCombination used to group the reductions relative to the various Jacobians (in general, more than one in transient problems) in a smart way. We consider, for example, the ODE

\[\tfrac{du}{dt} - \nu \Delta u = f \ \ \text{in} \ \ Ω \times [0,T]\]

subject to initial/boundary conditions. Upon applying a FE discretization in space, and a θ-method in time, one gets the space-time system

\[A_{\theta} u_{\theta} = f_{\theta}\]

where

\[A_{\theta} = \begin{bmatrix} A_1 + M / (\theta \Delta t) & & & & & \\ - M / (\theta \Delta t) & A_2 + M / (\theta \Delta t) & & & & \\ & - M / (\theta \Delta t) & A_3 + M / (\theta \Delta t) & & & \\ & & \ddots & \ddots & & \\ & & & & - M / (\theta \Delta t) & A_n + M / (\theta \Delta t) \end{bmatrix};\]

\[u_{\theta} = \begin{bmatrix} (1-\theta)u_0 + \theta u_1 \\ \vdots \\ (1-\theta)u_{n-1} + \theta u_n \end{bmatrix};\]

\[f_{\theta} = \begin{bmatrix} f_1 \\ \vdots \\ f_n \end{bmatrix};\]

\[A_k = A(t_{k-1} + \theta \Delta t);\]

\[f_k = f(t_{k-1} + \theta \Delta t).\]

Note: instead of multiplying $A_{\theta}$ by $u_{\theta}$, we multiply $\tilde{A}_{\theta}$ by $u$, where

\[\tilde{A}_{\theta} = tridiag((1-\theta)A_{k-1} - M / \Delta t, \theta A_k + M / \Delta t, 0).\]

We now denote with $\Phi$ and $\Psi$ the spatial and temporal basis obtained by reducing the snapshots associated to the state variable $u$. The Galerkin projection of the space-time system is equal to $\hat{A}_{\theta}\hat{u} = \hat{f}_{\theta}$, where $\hat{u}$ is the unknown, and

\[\begin{align*} \hat{A}_{\theta} &= \sum\limits_{k=1}^{n-1} ( (1-θ) \Phi^T A_k \Phi - \Phi^T M \Phi / \Delta t) \otimes \Psi[k-1,:]^T \Psi[k,:] + \sum\limits_{k=1}^n (\theta \Phi^T A_k \Phi + \Phi^T M \Phi / \Delta t) \otimes \Psi[k,:]^T \Psi[k,:] \\ &= \theta A_{backwards} + (1-\theta)A_{forwards} + (M_{backwards} + M_{forwards}) / \Delta t \\ \hat{f}_{\theta} &= \sum\limits_{k=1}^n \Phi^T f_k \otimes \Psi[k,:] \end{align*}\]

We notice that the expression of $\hat{A}_{\theta}$ can be written in a more general form as

\[\hat{A}_{\theta} = combine_A(A_{backwards},A_{forwards}) + combine_M(M_{backwards},M_{forwards}),\]

where $combine_A$ and $combine_M$ are two function specific to $A$ and $M$:

\[\begin{align*} combine_A(x,y) &= \theta y + (1-\theta)y \\ combine_M(x,y) &= (x - y) / \Delta t \end{align*}\]

The same can be said of any time marching scheme. This is the meaning of combine.

source
GridapROMs.RBTransient.KroneckerProjectionType
struct KroneckerProjection <: TransientProjection
  projection_space::Projection
  projection_time::Projection
end

Projection 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

source
GridapROMs.RBTransient.KroneckerReductionType
struct KroneckerReduction{A,B} <: HighDimReduction{A,B}
  reductions::AbstractVector{<:Reduction}
end

Wrapper for reduction methods in high order problems, such as transient ones. The reduced subspaces are constructed as Kronecker product spaces

source
GridapROMs.RBTransient.SequentialReductionType
struct SequentialReduction{A,B} <: HighDimReduction{A,B}
  reduction::Reduction{A,B}
end

Wrapper for sequential reduction methods in high-order problems, e.g. TT-SVD in transient applications

source
GridapROMs.RBTransient.SteadyReductionType
struct SteadyReduction{A,B} <: HighDimReduction{A,B}
  reduction::Reduction{A,B}
end

Wrapper 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.

source
GridapROMs.RBTransient.TransientIntegrationDomainType
struct TransientIntegrationDomain{A<:TransientIntegrationDomainStyle,Ti<:Integer} <: IntegrationDomain
  domain_style::A
  domain_space::IntegrationDomain
  indices_time::Vector{Ti}
end

Integration domain for a projection operator in a transient problem

source
GridapROMs.RBTransient.time_enrichmentMethod
time_enrichment(red::SupremizerReduction,a_primal::Projection,basis_dual;kwargs...) -> AbstractMatrix

Temporal supremizer enrichment. (Approximate) Procedure:

  1. for every b_dual ∈ Col(basis_dual)
  2. compute Φ_primal_dual = get_basis(a_primal)'*get_basis(b_dual)
  3. compute v = kernel(Φ_primal_dual)
  4. compute v′ = orth_complement(v,a_primal)
  5. enrich a_primal = [a_primal,v′]
source