Model Reduction

Balanced Truncation

resolvent4py.model_reduction.balanced_truncation.assemble_reduced_order_tensors(L: LinearOperator, B: BV, C: BV, Phi: BV, Psi: BV) Tuple[ndarray, ndarray, ndarray][source]

Assemble the reduced-order tensors \(A_r = \Psi^\intercal L \Phi\), \(B_r = \Psi^\intercal B\) and \(C_r = \Phi^\intercal C\).

Parameters:
  • L (LinearOperator) – linear operator governing the dynamics of the full-order system

  • B (SLEPc.BV) – input matrix of the full-order system

  • C (SLEPc.BV) – (complex-conj. transpose) of the output matrix of the full-order system

Returns:

Tuple with \(A_r\), \(B_r\) and \(C_r\)

Return type:

Tuple[numpy.ndarray, numpy.ndarray, numpy.ndarray]

resolvent4py.model_reduction.balanced_truncation.compute_balanced_projection(X: BV, Y: BV, r: int) Tuple[BV, BV, ndarray][source]

Given the output \((X, Y)\) of compute_gramian_factors(), compute \(\Phi\) and \(\Psi\) (each of dimension \(N\times r\)). Given the singular value decomposition \(Y^*X = U\Sigma V^*\), these are given by

\[\Phi = X V \Sigma^{-1/2},\quad \Psi = Y U \Sigma^{-1/2}.\]
Parameters:
  • X (PETSc.Mat) – reachability Gramian factor

  • Y (PETSc. Mat) – observability Gramian factor

  • r (int) – number of columns of \(\Phi\) and \(\Psi\)

Returns:

tuple \((\Phi, \Psi, \Sigma)\)

Return type:

Tuple[SLEPc.BV, SLEPc.BV, numpy.ndarray]

resolvent4py.model_reduction.balanced_truncation.compute_gramian_factors(L_generators: Callable, frequencies: ndarray, weights: ndarray, B: BV, C: BV) Tuple[Mat, Mat][source]

Compute the Gramian factors as outlined in section 2 of [Dergham2011]. In particular, we approximate the Reachability and Observability Gramians as follows

\[\begin{split}\begin{align} G_R &= \frac{1}{2\pi}\int_{-\infty}^{\infty}R(\omega)BB^*R(\omega)^*\ \,d\omega \approx \sum_j w_j X(\omega_j)X(\omega_j)^* = XX^* \\ G_O &= \frac{1}{2\pi}\int_{-\infty}^{\infty}R(\omega)^*CC^*R(\omega)\,\ d\omega \approx \sum_j w_j Y(\omega_j)Y(\omega_j)^* = YY^* \end{align}\end{split}\]

where \(\omega_j\in\Omega\) are quadrature points and \(w_j\) are corresponding quadrature weights. Here, \(X(\omega_j) = R(\omega_j)B\), and \(R(\omega) = \left(i\omega I - A\right)^{-1}\). When the linear operator \(A\) is real valued, the Gramians may be further manipulated into

\[G_R = \frac{1}{\pi}\int_{0}^{\infty}R(\omega)BB^*R(\omega)^*\ \,d\omega \approx \sum_j \delta_j w_j\left(X_\mathrm{r}(\omega_i)\ X_\mathrm{r}(\omega)^* + X_{\mathrm{i}}(\omega_j)\ X_{\mathrm{i}}(\omega_j)^*\right),\]

where the subscripts “r” and “i” denote the real and imaginary parts, and \(\delta_j = 1\) if \(\omega_j = 0\) and \(\delta_j = 1/2\) otherwise.

Parameters:
  • L_generators (Tuple[[Callable, Callable, Tuple[Callable, Callable, ...]]]) – tuple of functions that define the linear operator \(R(\omega)^{-1}\), the action of \(R(\omega)\) on matrices, and destructors to avoid memory leaks. (See examples.)

  • frequencies (numpy.ndarray) – quadrature points \(\omega_j\)

  • weights (numpy.ndarray) – quadrature weights \(w_j\) (see description above to see how the definition of \(w_j\) changes depending on whether the system is real-valued).

  • B (SLEPc.BV) – input matrix

  • C (SLEPc.BV) – output matrix

Returns:

tuple \((X, Y)\) (see definitions above)

Return type:

Tuple[PETSc.Mat, PETSc.Mat]

References

[Dergham2011]

Dergham et al., Model reduction for fluids using frequential snapshots, Physics of Fluids, 2011