pymor.reductors package

Submodules

basic module


class pymor.reductors.basic.GenericPGReductor(fom, W, V, bases_are_biorthonormal, vector_ranged_operators=('initial_data', ), product=None)[source]

Bases: pymor.core.interfaces.BasicInterface

Generic Petrov-Galerkin reductor.

Replaces each Operator of the given Model with the projection onto the span of the given projection matrices.

Parameters

fom
The Model which is to be reduced.
W
VectorArray containing the left projection matrix.
V
VectorArray containing the right projection matrix.
bases_are_biorthonormal
Indicate whether or not V and W are biorthonormal w.r.t. product.
vector_ranged_operators
List of keys in fom.operators for which the corresponding Operator should be biorthogonally projected (i.e. operators which map to vectors in contrast to bilinear forms which map to functionals).
product
Inner product for the projection of the Operators given by vector_ranged_operators.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce()[source]

Perform the Petrov-Galerkin projection.

Returns

The reduced Model.


class pymor.reductors.basic.GenericRBReductor(fom, RB=None, basis_is_orthonormal=None, vector_ranged_operators=('initial_data', ), product=None)[source]

Bases: pymor.core.interfaces.BasicInterface

Generic reduced basis reductor.

Replaces each Operator of the given Model with the Galerkin projection onto the span of the given reduced basis.

Parameters

fom
The Model which is to be reduced.
RB
VectorArray containing the reduced basis on which to project.
basis_is_orthonormal
If RB is specified, indicate whether or not the basis is orthonormal w.r.t. product.
vector_ranged_operators
List of keys in fom.operators for which the corresponding Operator should be orthogonally projected (i.e. operators which map to vectors in contrast to bilinear forms which map to functionals).
product
Inner product for the orthonormalization of RB and the projection of the Operators given by vector_ranged_operators.
extend_basis(U, method='gram_schmidt', pod_modes=1, pod_orthonormalize=True, orthonormal=None, copy_U=True)[source]

Extend basis by new vectors.

Parameters

U
VectorArray containing the new basis vectors.
method

Basis extension method to use. The following methods are available:

trivial:Vectors in U are appended to the basis. Duplicate vectors in the sense of almost_equal are removed.
gram_schmidt:New basis vectors are orthonormalized w.r.t. to the old basis using the gram_schmidt algorithm.
pod:Append the first POD modes of the defects of the projections of the vectors in U onto the existing basis (e.g. for use in POD-Greedy algorithm).

Warning

In case of the 'gram_schmidt' and 'pod' extension methods, the existing reduced basis is assumed to be orthonormal w.r.t. the given inner product.

pod_modes
In case method == 'pod', the number of POD modes that shall be appended to the basis.
pod_orthonormalize
If True and method == 'pod', re-orthonormalize the new basis vectors obtained by the POD in order to improve numerical accuracy.
orthonormal
If method == 'trivial', set this to True to indicate that the basis will remain orthonormal after extending.
copy_U
If copy_U is False, the new basis vectors might be removed from U.

Raises

ExtensionError
Raised when the selected extension method does not yield a basis of increased dimension.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(dim=None)[source]

Perform the reduced basis projection.

Parameters

dim
If specified, the desired reduced state dimension. Must not be larger than the current reduced basis dimension.

Returns

The reduced Model.


pymor.reductors.basic.extend_basis(basis, U, product=None, method='gram_schmidt', pod_modes=1, pod_orthonormalize=True, copy_U=True)[source]

bt module


class pymor.reductors.bt.BRBTReductor(fom, gamma, solver_options=None)[source]

Bases: pymor.reductors.bt.GenericBTReductor

Bounded Real (BR) Balanced Truncation reductor.

See [A05] (Section 7.5.3) and [OJ88].

Parameters

fom
The system which is to be reduced.
gamma
Upper bound for the \mathcal{H}_\infty-norm.
solver_options
The solver options to use to solve the positive Riccati equations.
error_bounds()[source]

Returns error bounds for all possible reduced orders.

gramians()[source]

Return low-rank Cholesky factors of Gramians.


class pymor.reductors.bt.BTReductor(fom)[source]

Bases: pymor.reductors.bt.GenericBTReductor

Standard (Lyapunov) Balanced Truncation reductor.

See Section 7.3 in [A05].

Parameters

fom
The system which is to be reduced.
error_bounds()[source]

Returns error bounds for all possible reduced orders.

gramians()[source]

Return low-rank Cholesky factors of Gramians.


class pymor.reductors.bt.GenericBTReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Generic Balanced Truncation reductor.

Parameters

fom
The system which is to be reduced.
error_bounds()[source]

Returns error bounds for all possible reduced orders.

gramians()[source]

Return low-rank Cholesky factors of Gramians.

reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(r=None, tol=None, projection='bfsr')[source]

Generic Balanced Truncation.

Parameters

r
Order of the reduced model if tol is None, maximum order if tol is specified.
tol
Tolerance for the error bound if r is None.
projection

Projection method used:

  • 'sr': square root method
  • 'bfsr': balancing-free square root method (default, since it avoids scaling by singular values and orthogonalizes the projection matrices, which might make it more accurate than the square root method)
  • 'biorth': like the balancing-free square root method, except it biorthogonalizes the projection matrices (using gram_schmidt_biorth)

Returns

rom
Reduced system.
sv_U_V()[source]

Return singular values and vectors.


class pymor.reductors.bt.LQGBTReductor(fom, solver_options=None)[source]

Bases: pymor.reductors.bt.GenericBTReductor

Linear Quadratic Gaussian (LQG) Balanced Truncation reductor.

See Section 3 in [MG91].

Parameters

fom
The system which is to be reduced.
solver_options
The solver options to use to solve the Riccati equations.
error_bounds()[source]

Returns error bounds for all possible reduced orders.

gramians()[source]

Return low-rank Cholesky factors of Gramians.

coercive module


class pymor.reductors.coercive.CoerciveRBEstimator(residual, residual_range_dims, coercivity_estimator)[source]

Bases: pymor.core.interfaces.ImmutableInterface

Instantiated by CoerciveRBReductor.

Not to be used directly.


class pymor.reductors.coercive.CoerciveRBReductor(fom, RB=None, basis_is_orthonormal=None, vector_ranged_operators=('initial_data', ), product=None, coercivity_estimator=None)[source]

Bases: pymor.reductors.basic.GenericRBReductor

Reduced Basis reductor for StationaryModels with coercive linear operator.

The only addition to GenericRBReductor is an error estimator which evaluates the dual norm of the residual with respect to a given inner product. For the reduction of the residual we use ResidualReductor for improved numerical stability [BEOR14].

Parameters

fom
The Model which is to be reduced.
RB
VectorArray containing the reduced basis on which to project.
basis_is_orthonormal
If RB is specified, indicate whether or not the basis is orthonormal w.r.t. product.
vector_ranged_operators
List of keys in fom.operators for which the corresponding Operator should be orthogonally projected (i.e. operators which map to vectors in contrast to bilinear forms which map to functionals).
product
Inner product for the orthonormalization of RB, the projection of the Operators given by vector_ranged_operators and for the computation of Riesz representatives of the residual. If None, the Euclidean product is used.
coercivity_estimator
None or a ParameterFunctional returning a lower bound for the coercivity constant of the given problem. Note that the computed error estimate is only guaranteed to be an upper bound for the error when an appropriate coercivity estimate is specified.

class pymor.reductors.coercive.SimpleCoerciveRBEstimator(estimator_matrix, coercivity_estimator)[source]

Bases: pymor.core.interfaces.ImmutableInterface

Instantiated by SimpleCoerciveRBReductor.

Not to be used directly.


class pymor.reductors.coercive.SimpleCoerciveRBReductor(fom, RB=None, basis_is_orthonormal=None, vector_ranged_operators=('initial_data', ), product=None, coercivity_estimator=None)[source]

Bases: pymor.reductors.basic.GenericRBReductor

Reductor for linear StationaryModels with affinely decomposed operator and rhs.

Note

The reductor CoerciveRBReductor can be used for arbitrary coercive StationaryModels and offers an improved error estimator with better numerical stability.

The only addition is to GenericRBReductor is an error estimator, which evaluates the norm of the residual with respect to a given inner product.

Parameters

fom
The Model which is to be reduced.
RB
VectorArray containing the reduced basis on which to project.
basis_is_orthonormal
If RB is specified, indicate whether or not the basis is orthonormal w.r.t. product.
vector_ranged_operators
List of keys in fom.operators for which the corresponding Operator should be orthogonally projected (i.e. operators which map to vectors in contrast to bilinear forms which map to functionals).
product
Inner product for the orthonormalization of RB, the projection of the Operators given by vector_ranged_operators and for the computation of Riesz representatives of the residual. If None, the Euclidean product is used.
coercivity_estimator
None or a ParameterFunctional returning a lower bound for the coercivity constant of the given problem. Note that the computed error estimate is only guaranteed to be an upper bound for the error when an appropriate coercivity estimate is specified.

h2 module


class pymor.reductors.h2.IRKAReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Iterative Rational Krylov Algorithm reductor.

Parameters

fom
LTIModel.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(r, sigma=None, b=None, c=None, rom0=None, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, projection='orth', use_arnoldi=False, conv_crit='sigma', compute_errors=False)[source]

Reduce using IRKA.

See [GAB08] (Algorithm 4.1) and [ABG10] (Algorithm 1).

Parameters

r
Order of the reduced order model.
sigma

Initial interpolation points (closed under conjugation).

If None, interpolation points are log-spaced between 0.1 and 10. If sigma is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a one-dimensional array-like of length r.

sigma and rom0 cannot both be not None.

b

Initial right tangential directions.

If None, if is chosen as all ones. If b is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a VectorArray of length r from fom.B.source.

b and rom0 cannot both be not None.

c

Initial left tangential directions.

If None, if is chosen as all ones. If c is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a VectorArray of length r from fom.C.range.

c and rom0 cannot both be not None.

rom0

Initial reduced order model.

If None, then sigma, b, and c are used. Otherwise, it needs to be an LTIModel of order r and it is used to construct sigma, b, and c.

tol
Tolerance for the convergence criterion.
maxit
Maximum number of iterations.
num_prev
Number of previous iterations to compare the current iteration to. Larger number can avoid occasional cyclic behavior of IRKA.
force_sigma_in_rhp
If False, new interpolation are reflections of the current reduced order model’s poles. Otherwise, only poles in the left half-plane are reflected.
projection

Projection method:

  • 'orth': projection matrices are orthogonalized with respect to the Euclidean inner product
  • 'biorth': projection matrices are biorthogolized with respect to the E product
use_arnoldi
Should the Arnoldi process be used for rational interpolation. Available only for SISO systems. Otherwise, it is ignored.
conv_crit

Convergence criterion:

  • 'sigma': relative change in interpolation points
  • 'h2': relative \mathcal{H}_2 distance of reduced-order models
compute_errors

Should the relative \mathcal{H}_2-errors of intermediate reduced order models be computed.

Warning

Computing \mathcal{H}_2-errors is expensive. Use this option only if necessary.

Returns

rom
Reduced LTIModel model.

class pymor.reductors.h2.OneSidedIRKAReductor(fom, version)[source]

Bases: pymor.core.interfaces.BasicInterface

One-Sided Iterative Rational Krylov Algorithm reductor.

Parameters

fom
LTIModel.
version

Version of the one-sided IRKA:

  • 'V': Galerkin projection using the input Krylov subspace,
  • 'W': Galerkin projection using the output Krylov subspace.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(r, sigma=None, b=None, c=None, rd0=None, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, projection='orth', conv_crit='sigma', compute_errors=False)[source]

Reduce using one-sided IRKA.

Parameters

r
Order of the reduced order model.
sigma

Initial interpolation points (closed under conjugation).

If None, interpolation points are log-spaced between 0.1 and 10. If sigma is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a one-dimensional array-like of length r.

sigma and rd0 cannot both be not None.

b

Initial right tangential directions.

If None, if is chosen as all ones. If b is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a VectorArray of length r from fom.B.source.

b and rd0 cannot both be not None.

c

Initial left tangential directions.

If None, if is chosen as all ones. If c is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a VectorArray of length r from fom.C.range.

c and rd0 cannot both be not None.

rd0

Initial reduced order model.

If None, then sigma, b, and c are used. Otherwise, it needs to be an LTIModel of order r and it is used to construct sigma, b, and c.

tol
Tolerance for the largest change in interpolation points.
maxit
Maximum number of iterations.
num_prev
Number of previous iterations to compare the current iteration to. A larger number can avoid occasional cyclic behavior.
force_sigma_in_rhp
If ‘False`, new interpolation are reflections of reduced order model’s poles. Otherwise, they are always in the right half-plane.
projection

Projection method:

  • 'orth': projection matrix is orthogonalized with respect to the Euclidean inner product,
  • 'Eorth': projection matrix is orthogonalized with respect to the E product.
conv_crit

Convergence criterion:

  • 'sigma': relative change in interpolation points,
  • 'h2': relative \mathcal{H}_2 distance of reduced order models.
compute_errors

Should the relative \mathcal{H}_2-errors of intermediate reduced order models be computed.

Warning

Computing \mathcal{H}_2-errors is expensive. Use this option only if necessary.

Returns

rom
Reduced LTIModel model.

class pymor.reductors.h2.TF_IRKAReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Realization-independent IRKA reductor.

See [BG12].

Parameters

fom
Model with eval_tf and eval_dtf methods.
reduce(r, sigma=None, b=None, c=None, rom0=None, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, conv_crit='sigma')[source]

Reduce using TF-IRKA.

Parameters

r
Order of the reduced order model.
sigma

Initial interpolation points (closed under conjugation).

If None, interpolation points are log-spaced between 0.1 and 10. If sigma is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a one-dimensional array-like of length r.

sigma and rom0 cannot both be not None.

b

Initial right tangential directions.

If None, if is chosen as all ones. If b is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a NumPy array of shape (m, r).

b and rom0 cannot both be not None.

c

Initial left tangential directions.

If None, if is chosen as all ones. If c is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a NumPy array of shape (p, r).

c and rom0 cannot both be not None.

rom0

Initial reduced order model.

If None, then sigma, b, and c are used. Otherwise, it needs to be an LTIModel of order r and it is used to construct sigma, b, and c.

tol
Tolerance for the convergence criterion.
maxit
Maximum number of iterations.
num_prev
Number of previous iterations to compare the current iteration to. Larger number can avoid occasional cyclic behavior of TF-IRKA.
force_sigma_in_rhp
If False, new interpolation are reflections of the current reduced order model’s poles. Otherwise, only the poles in the left half-plane are reflected.
conv_crit

Convergence criterion:

  • 'sigma': relative change in interpolation points
  • 'h2': relative \mathcal{H}_2 distance of reduced-order models

Returns

rom
Reduced LTIModel model.

class pymor.reductors.h2.TSIAReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Two-Sided Iteration Algorithm reductor.

Parameters

fom
LTIModel.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(rom0, tol=0.0001, maxit=100, num_prev=1, projection='orth', conv_crit='sigma', compute_errors=False)[source]

Reduce using TSIA.

See [XZ11] (Algorithm 1) and [BKS11].

In exact arithmetic, TSIA is equivalent to IRKA (under some assumptions on the poles of the reduced model). The main difference in implementation is that TSIA computes the Schur decomposition of the reduced matrices, while IRKA computes the eigenvalue decomposition. Therefore, TSIA might behave better for non-normal reduced matrices.

Parameters

rom0
Initial reduced order model.
tol
Tolerance for the convergence criterion.
maxit
Maximum number of iterations.
num_prev
Number of previous iterations to compare the current iteration to. Larger number can avoid occasional cyclic behavior of TSIA.
projection

Projection method:

  • 'orth': projection matrices are orthogonalized with respect to the Euclidean inner product
  • 'biorth': projection matrices are biorthogolized with respect to the E product
conv_crit

Convergence criterion:

  • 'sigma': relative change in interpolation points
  • 'h2': relative \mathcal{H}_2 distance of reduced-order models
compute_errors

Should the relative \mathcal{H}_2-errors of intermediate reduced order models be computed.

Warning

Computing \mathcal{H}_2-errors is expensive. Use this option only if necessary.

Returns

rom
Reduced LTIModel.

pymor.reductors.h2._convergence_criterion(data, conv_crit)[source]

Compute the convergence criterion for given data.


pymor.reductors.h2._poles_and_tangential_directions(rom)[source]

Compute the poles and tangential directions of a reduced order model.

interpolation module


class pymor.reductors.interpolation.DelayBHIReductor(fom)[source]

Bases: pymor.reductors.interpolation.GenericBHIReductor

Bitangential Hermite interpolation for delay systems.

Parameters

fom
LinearDelayModel.

class pymor.reductors.interpolation.GenericBHIReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Generic bitangential Hermite interpolation reductor.

This is a generic reductor for reducing any linear InputStateOutputModel with the transfer function which can be written in the generalized coprime factorization \mathcal{C}(s) \mathcal{K}(s)^{-1}
\mathcal{B}(s) as in [BG09]. The interpolation here is limited to only up to the first derivative. Hence, interpolation points are assumed to be pairwise distinct.

Parameters

fom
Model.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(sigma, b, c, projection='orth')[source]

Bitangential Hermite interpolation.

Parameters

sigma
Interpolation points (closed under conjugation), list of length r.
b
Right tangential directions, VectorArray of length r from self.fom.input_space.
c
Left tangential directions, VectorArray of length r from self.fom.output_space.
projection

Projection method:

  • 'orth': projection matrices are orthogonalized with respect to the Euclidean inner product
  • 'biorth': projection matrices are biorthogolized with respect to the E product

Returns

rom
Reduced model.

class pymor.reductors.interpolation.LTI_BHIReductor(fom)[source]

Bases: pymor.reductors.interpolation.GenericBHIReductor

Bitangential Hermite interpolation for LTIModels.

Parameters

fom
LTIModel.
reduce(sigma, b, c, projection='orth', use_arnoldi=False)[source]

Bitangential Hermite interpolation.

Parameters

sigma
Interpolation points (closed under conjugation), list of length r.
b
Right tangential directions, VectorArray of length r from self.fom.input_space.
c
Left tangential directions, VectorArray of length r from self.fom.output_space.
projection

Projection method:

  • 'orth': projection matrices are orthogonalized with respect to the Euclidean inner product
  • 'biorth': projection matrices are biorthogolized with respect to the E product
use_arnoldi
Should the Arnoldi process be used for rational interpolation. Available only for SISO systems. Otherwise, it is ignored.

Returns

rom
Reduced model.
reduce_arnoldi(sigma, b, c)[source]

Bitangential Hermite interpolation for SISO LTIModels.

Parameters

sigma
Interpolation points (closed under conjugation), list of length r.
b
Right tangential directions, VectorArray of length r from self.fom.B.source.
c
Left tangential directions, VectorArray of length r from self.fom.C.range.

Returns

rom
Reduced LTIModel model.

class pymor.reductors.interpolation.SO_BHIReductor(fom)[source]

Bases: pymor.reductors.interpolation.GenericBHIReductor

Bitangential Hermite interpolation for second-order systems.

Parameters

fom
SecondOrderModel.

class pymor.reductors.interpolation.TFInterpReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Loewner bitangential Hermite interpolation reductor.

See [BG12].

Parameters

fom
Model with eval_tf and eval_dtf methods.
reduce(sigma, b, c)[source]

Realization-independent tangential Hermite interpolation.

Parameters

sigma
Interpolation points (closed under conjugation), list of length r.
b
Right tangential directions, NumPy array of shape (fom.input_dim, r).
c
Left tangential directions, NumPy array of shape (fom.output_dim, r).

Returns

lti
LTIModel interpolating the transfer function of fom.

parabolic module


class pymor.reductors.parabolic.ParabolicRBEstimator(residual, residual_range_dims, initial_residual, initial_residual_range_dims, coercivity_estimator)[source]

Bases: pymor.core.interfaces.ImmutableInterface

Instantiated by ParabolicRBReductor.

Not to be used directly.


class pymor.reductors.parabolic.ParabolicRBReductor(fom, RB=None, basis_is_orthonormal=None, product=None, coercivity_estimator=None)[source]

Bases: pymor.reductors.basic.GenericRBReductor

Reduced Basis Reductor for parabolic equations.

This reductor uses GenericRBReductor for the actual RB-projection. The only addition is the assembly of an error estimator which bounds the discrete l2-in time / energy-in space error similar to [GP05], [HO08] as follows:

\left[ C_a^{-1}(\mu)\|e_N(\mu)\|^2 + \sum_{n=1}^{N} \Delta t\|e_n(\mu)\|^2_e \right]^{1/2}
    \leq \left[ C_a^{-1}(\mu)\Delta t \sum_{n=1}^{N}\|\mathcal{R}^n(u_n(\mu), \mu)\|^2_{e,-1}
                + C_a^{-1}(\mu)\|e_0\|^2 \right]^{1/2}

Here, \|\cdot\| denotes the norm induced by the problem’s mass matrix (e.g. the L^2-norm) and \|\cdot\|_e is an arbitrary energy norm w.r.t. which the space operator A(\mu) is coercive, and C_a(\mu) is a lower bound for its coercivity constant. Finally, \mathcal{R}^n denotes the implicit Euler timestepping residual for the (fixed) time step size \Delta t,

\mathcal{R}^n(u_n(\mu), \mu) :=
    f - M \frac{u_{n}(\mu) - u_{n-1}(\mu)}{\Delta t} - A(u_n(\mu), \mu),

where M denotes the mass operator and f the source term. The dual norm of the residual is computed using the numerically stable projection from [BEOR14].

Parameters

fom
The InstationaryModel which is to be reduced.
RB
VectorArray containing the reduced basis on which to project.
basis_is_orthonormal
Indicate whether or not the basis is orthonormal w.r.t. product.
product
The energy inner product Operator w.r.t. which the reduction error is estimated and RB is orthonormalized.
coercivity_estimator
None or a ParameterFunctional returning a lower bound C_a(\mu) for the coercivity constant of fom.operator w.r.t. product.

residual module


class pymor.reductors.residual.ImplicitEulerResidualOperator(operator, mass, rhs, dt, name=None)[source]

Bases: pymor.operators.basic.OperatorBase

Instantiated by ImplicitEulerResidualReductor.

apply(U, U_old, mu=None)[source]

Apply the operator to a VectorArray.

Parameters

U
VectorArray of vectors to which the operator is applied.
mu
The Parameter for which to evaluate the operator.

Returns

VectorArray of the operator evaluations.


class pymor.reductors.residual.ImplicitEulerResidualReductor(RB, operator, mass, dt, rhs=None, product=None)[source]

Bases: pymor.core.interfaces.BasicInterface

Reduced basis residual reductor with mass operator for implicit Euler timestepping.

Given an operator, mass and a functional, the concatenation of residual operator with the Riesz isomorphism is given by:

riesz_residual.apply(U, U_old, mu)
    == product.apply_inverse(operator.apply(U, mu) + 1/dt*mass.apply(U, mu) - 1/dt*mass.apply(U_old, mu)
       - rhs.as_vector(mu))

This reductor determines a low-dimensional subspace of the image of a reduced basis space under riesz_residual using estimate_image_hierarchical, computes an orthonormal basis residual_range of this range space and then returns the Petrov-Galerkin projection

projected_riesz_residual
    == riesz_residual.projected(range_basis=residual_range, source_basis=RB)

of the riesz_residual operator. Given reduced basis coefficient vectors u and u_old, the dual norm of the residual can then be computed as

projected_riesz_residual.apply(u, u_old, mu).l2_norm()

Moreover, a reconstruct method is provided such that

residual_reductor.reconstruct(projected_riesz_residual.apply(u, u_old, mu))
    == riesz_residual.apply(RB.lincomb(u), RB.lincomb(u_old), mu)

Parameters

operator
See definition of riesz_residual.
mass
The mass operator. See definition of riesz_residual.
dt
The time step size. See definition of riesz_residual.
rhs
See definition of riesz_residual. If None, zero right-hand side is assumed.
RB
VectorArray containing a basis of the reduced space onto which to project.
product
Inner product Operator w.r.t. which to compute the Riesz representatives.
reconstruct(u)[source]

Reconstruct high-dimensional residual vector from reduced vector u.


class pymor.reductors.residual.NonProjectedImplicitEulerResidualOperator(operator, mass, rhs, dt, product)[source]

Bases: pymor.reductors.residual.ImplicitEulerResidualOperator

Instantiated by ImplicitEulerResidualReductor.

Not to be used directly.

apply(U, U_old, mu=None)[source]

Apply the operator to a VectorArray.

Parameters

U
VectorArray of vectors to which the operator is applied.
mu
The Parameter for which to evaluate the operator.

Returns

VectorArray of the operator evaluations.


class pymor.reductors.residual.NonProjectedResidualOperator(operator, rhs, riesz_representatives, product)[source]

Bases: pymor.reductors.residual.ResidualOperator

Instantiated by ResidualReductor.

Not to be used directly.

apply(U, mu=None)[source]

Apply the operator to a VectorArray.

Parameters

U
VectorArray of vectors to which the operator is applied.
mu
The Parameter for which to evaluate the operator.

Returns

VectorArray of the operator evaluations.


class pymor.reductors.residual.ResidualOperator(operator, rhs, name=None)[source]

Bases: pymor.operators.basic.OperatorBase

Instantiated by ResidualReductor.

apply(U, mu=None)[source]

Apply the operator to a VectorArray.

Parameters

U
VectorArray of vectors to which the operator is applied.
mu
The Parameter for which to evaluate the operator.

Returns

VectorArray of the operator evaluations.


class pymor.reductors.residual.ResidualReductor(RB, operator, rhs=None, product=None, riesz_representatives=False)[source]

Bases: pymor.core.interfaces.BasicInterface

Generic reduced basis residual reductor.

Given an operator and a right-hand side, the residual is given by:

residual.apply(U, mu) == operator.apply(U, mu) - rhs.as_range_array(mu)

When operator maps to functionals instead of vectors, we are interested in the Riesz representative of the residual:

residual.apply(U, mu)
    == product.apply_inverse(operator.apply(U, mu) - rhs.as_range_array(mu))

Given a basis RB of a subspace of the source space of operator, this reductor uses estimate_image_hierarchical to determine a low-dimensional subspace containing the image of the subspace under residual (resp. riesz_residual), computes an orthonormal basis residual_range for this range space and then returns the Petrov-Galerkin projection

projected_residual
    == project(residual, range_basis=residual_range, source_basis=RB)

of the residual operator. Given a reduced basis coefficient vector u, w.r.t. RB, the (dual) norm of the residual can then be computed as

projected_residual.apply(u, mu).l2_norm()

Moreover, a reconstruct method is provided such that

residual_reductor.reconstruct(projected_residual.apply(u, mu))
    == residual.apply(RB.lincomb(u), mu)

Parameters

RB
VectorArray containing a basis of the reduced space onto which to project.
operator
See definition of residual.
rhs
See definition of residual. If None, zero right-hand side is assumed.
product
Inner product Operator w.r.t. which to orthonormalize and w.r.t. which to compute the Riesz representatives in case operator maps to functionals.
riesz_representatives
If True compute the Riesz representative of the residual.
reconstruct(u)[source]

Reconstruct high-dimensional residual vector from reduced vector u.

sobt module


class pymor.reductors.sobt.GenericSOBTpvReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Generic Second-Order Balanced Truncation position/velocity reductor.

See [RS08].

Parameters

fom
The system which is to be reduced.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(r, projection='bfsr')[source]

Reduce using GenericSOBTpv.

Parameters

r
Order of the reduced model.
projection

Projection method used:

  • 'sr': square root method
  • 'bfsr': balancing-free square root method (default, since it avoids scaling by singular values and orthogonalizes the projection matrices, which might make it more accurate than the square root method)
  • 'biorth': like the balancing-free square root method, except it biorthogonalizes the projection matrices

Returns

rom
Reduced system.

class pymor.reductors.sobt.SOBTReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Second-Order Balanced Truncation reductor.

See [CLVV06].

Parameters

fom
The system which is to be reduced.
reduce(r, projection='bfsr')[source]

Reduce using SOBT.

Parameters

r
Order of the reduced model.
projection

Projection method used:

  • 'sr': square root method
  • 'bfsr': balancing-free square root method (default, since it avoids scaling by singular values and orthogonalizes the projection matrices, which might make it more accurate than the square root method)
  • 'biorth': like the balancing-free square root method, except it biorthogonalizes the projection matrices

Returns

rom
Reduced system.

class pymor.reductors.sobt.SOBTfvReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

Free-velocity Second-Order Balanced Truncation reductor.

See [MS96].

Parameters

fom
The system which is to be reduced.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(r, projection='bfsr')[source]

Reduce using SOBTfv.

Parameters

r
Order of the reduced model.
projection

Projection method used:

  • 'sr': square root method
  • 'bfsr': balancing-free square root method (default, since it avoids scaling by singular values and orthogonalizes the projection matrices, which might make it more accurate than the square root method)
  • 'biorth': like the balancing-free square root method, except it biorthogonalizes the projection matrices

Returns

rom
Reduced system.

class pymor.reductors.sobt.SOBTpReductor(fom)[source]

Bases: pymor.reductors.sobt.GenericSOBTpvReductor

Second-Order Balanced Truncation position reductor.

See [RS08].

Parameters

fom
The system which is to be reduced.

class pymor.reductors.sobt.SOBTpvReductor(fom)[source]

Bases: pymor.reductors.sobt.GenericSOBTpvReductor

Second-Order Balanced Truncation position-velocity reductor.

See [RS08].

Parameters

fom
The system which is to be reduced.

class pymor.reductors.sobt.SOBTvReductor(fom)[source]

Bases: pymor.reductors.sobt.GenericSOBTpvReductor

Second-Order Balanced Truncation velocity reductor.

See [RS08].

Parameters

fom
The system which is to be reduced.

class pymor.reductors.sobt.SOBTvpReductor(fom)[source]

Bases: pymor.reductors.sobt.GenericSOBTpvReductor

Second-Order Balanced Truncation velocity-position reductor.

See [RS08].

Parameters

fom
The system which is to be reduced.

sor_irka module


class pymor.reductors.sor_irka.SOR_IRKAReductor(fom)[source]

Bases: pymor.core.interfaces.BasicInterface

SOR-IRKA reductor.

Parameters

fom
SecondOrderModel.
reconstruct(u)[source]

Reconstruct high-dimensional vector from reduced vector u.

reduce(r, sigma=None, b=None, c=None, rom0=None, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, projection='orth', use_arnoldi=False, conv_crit='sigma', compute_errors=False, irka_options=None)[source]

Reduce using SOR-IRKA.

It uses IRKA as the intermediate reductor, to reduce from 2r to r poles. See Section 5.3.2 in [W12].

Parameters

r
Order of the reduced order model.
sigma

Initial interpolation points (closed under conjugation).

If None, interpolation points are log-spaced between 0.1 and 10. If sigma is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a one-dimensional array-like of length r.

sigma and rom0 cannot both be not None.

b

Initial right tangential directions.

If None, if is chosen as all ones. If b is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a VectorArray of length r from fom.B.source.

b and rom0 cannot both be not None.

c

Initial left tangential directions.

If None, if is chosen as all ones. If c is an int, it is used as a seed to generate it randomly. Otherwise, it needs to be a VectorArray of length r from fom.Cp.range.

c and rom0 cannot both be not None.

rom0

Initial reduced order model.

If None, then sigma, b, and c are used. Otherwise, it needs to be an LTIModel of order r and it is used to construct sigma, b, and c.

tol
Tolerance for the convergence criterion.
maxit
Maximum number of iterations.
num_prev
Number of previous iterations to compare the current iteration to. Larger number can avoid occasional cyclic behavior of IRKA.
force_sigma_in_rhp
If False, new interpolation are reflections of the current reduced order model’s poles. Otherwise, only the poles in the left half-plane are reflected.
projection

Projection method:

  • 'orth': projection matrices are orthogonalized with respect to the Euclidean inner product
  • 'biorth': projection matrices are biorthogolized with respect to the E product
conv_crit

Convergence criterion:

  • 'sigma': relative change in interpolation points
  • 'h2': relative \mathcal{H}_2 distance of reduced-order models
compute_errors

Should the relative \mathcal{H}_2-errors of intermediate reduced order models be computed.

Warning

Computing \mathcal{H}_2-errors is expensive. Use this option only if necessary.

irka_options
Dict of options for IRKAReductor.reduce.

Returns

rom
Reduced LTIModel model.