pymor.reductors package¶
Submodules¶
basic module¶
-
class
pymor.reductors.basic.DelayLTIPGReductor(fom, W, V, E_biorthonormal=False)[source]¶ Bases:
pymor.reductors.basic.ProjectionBasedReductorPetrov-Galerkin projection of an
LinearDelayModel.Parameters
Methods
build_rom,extend_basis,project_operators,project_operators_to_subbasis,reconstructassemble_estimator,assemble_estimator_for_subbasis,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.basic.InstationaryRBReductor(fom, RB=None, product=None, initial_data_product=None, product_is_mass=False, check_orthonormality=None, check_tol=None)[source]¶ Bases:
pymor.reductors.basic.ProjectionBasedReductorGalerkin projection of an
InstationaryModel.Parameters
- fom
The full order
Modelto reduce.- RB
The basis of the reduced space onto which to project. If
Nonean empty basis is used.- product
Inner product
Operatorw.r.t. whichRBis orthonormalized. IfNone, the the Euclidean inner product is used.- initial_data_product
Inner product
Operatorw.r.t. which theinitial_dataoffomis orthogonally projected. IfNone, the Euclidean inner product is used.- product_is_mass
If
True, no mass matrix for the reducedModelis assembled. Set toTrueifRBis orthonormal w.r.t. themassmatrix offom.- check_orthonormality
- check_tol
Methods
build_rom,project_operators,project_operators_to_subbasisassemble_estimator,assemble_estimator_for_subbasis,extend_basis,reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.basic.LTIPGReductor(fom, W, V, E_biorthonormal=False)[source]¶ Bases:
pymor.reductors.basic.ProjectionBasedReductorPetrov-Galerkin projection of an
LTIModel.Parameters
Methods
build_rom,extend_basis,project_operators,project_operators_to_subbasis,reconstructassemble_estimator,assemble_estimator_for_subbasis,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.basic.ProjectionBasedReductor(fom, bases, products={}, check_orthonormality=True, check_tol=0.001)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric projection based reductor.
Parameters
- fom
The full order
Modelto reduce.- bases
A dict of
VectorArraysof basis vectors.- products
A dict of inner product
Operatorsw.r.t. which the corresponding bases are orthonormalized. A value ofNonecorresponds to orthonormalization of the basis w.r.t. the Euclidean inner product.- check_orthonormality
If
True, check if bases which have a corresponding entry in theproductsdict are orthonormal w.r.t. the given inner product. After eachbasis extension, orthonormality is checked again.- check_tol
If
check_orthonormalityisTrue, the numerical tolerance with which the checks are performed.
Methods
assemble_estimator,assemble_estimator_for_subbasis,build_rom,extend_basis,project_operators,project_operators_to_subbasis,reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.basic.SOLTIPGReductor(fom, W, V, M_biorthonormal=False)[source]¶ Bases:
pymor.reductors.basic.ProjectionBasedReductorPetrov-Galerkin projection of an
SecondOrderModel.Parameters
Methods
build_rom,extend_basis,project_operators,project_operators_to_subbasis,reconstructassemble_estimator,assemble_estimator_for_subbasis,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.basic.StationaryRBReductor(fom, RB=None, product=None, check_orthonormality=None, check_tol=None)[source]¶ Bases:
pymor.reductors.basic.ProjectionBasedReductorGalerkin projection of a
StationaryModel.Parameters
Methods
build_rom,project_operators,project_operators_to_subbasisassemble_estimator,assemble_estimator_for_subbasis,extend_basis,reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
bt module¶
-
class
pymor.reductors.bt.BRBTReductor(fom, gamma=1, mu=None, solver_options=None)[source]¶ Bases:
pymor.reductors.bt.GenericBTReductorBounded Real (BR) Balanced Truncation reductor.
See [A05] (Section 7.5.3) and [OJ88].
Parameters
Attributes
-
class
pymor.reductors.bt.BTReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.bt.GenericBTReductorStandard (Lyapunov) Balanced Truncation reductor.
See Section 7.3 in [A05].
Attributes
-
class
pymor.reductors.bt.GenericBTReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric Balanced Truncation reductor.
Attributes
-
reduce(r=None, tol=None, projection='bfsr')[source]¶ Generic Balanced Truncation.
Parameters
- r
Order of the reduced model if
tolisNone, maximum order iftolis specified.- tol
Tolerance for the error bound if
risNone.- 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 (usinggram_schmidt_biorth)
Returns
- rom
Reduced-order model.
-
-
class
pymor.reductors.bt.LQGBTReductor(fom, mu=None, solver_options=None)[source]¶ Bases:
pymor.reductors.bt.GenericBTReductorLinear Quadratic Gaussian (LQG) Balanced Truncation reductor.
See Section 3 in [MG91].
Parameters
Attributes
coercive module¶
-
class
pymor.reductors.coercive.CoerciveRBEstimator(residual, residual_range_dims, coercivity_estimator)[source]¶ Bases:
pymor.core.interfaces.ImmutableInterfaceInstantiated by
CoerciveRBReductor.Not to be used directly.
Methods
estimate,restricted_to_subbasisdisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.coercive.CoerciveRBReductor(fom, RB=None, product=None, coercivity_estimator=None, check_orthonormality=None, check_tol=None)[source]¶ Bases:
pymor.reductors.basic.StationaryRBReductorReduced Basis reductor for
StationaryModelswith coercive linear operator.The only addition to
StationaryRBReductoris 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 useResidualReductorfor improved numerical stability [BEOR14].Parameters
- fom
The
Modelwhich is to be reduced.- RB
VectorArraycontaining the reduced basis on which to project.- product
Inner product for the orthonormalization of
RB, the projection of theOperatorsgiven byvector_ranged_operatorsand for the computation of Riesz representatives of the residual. IfNone, the Euclidean product is used.- coercivity_estimator
Noneor aParameterFunctionalreturning 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.
Methods
assemble_estimator,assemble_estimator_for_subbasisbuild_rom,project_operators,project_operators_to_subbasisextend_basis,reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.coercive.SimpleCoerciveRBEstimator(estimator_matrix, coercivity_estimator)[source]¶ Bases:
pymor.core.interfaces.ImmutableInterfaceInstantiated by
SimpleCoerciveRBReductor.Not to be used directly.
Methods
estimate,restricted_to_subbasisdisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.coercive.SimpleCoerciveRBReductor(fom, RB=None, product=None, coercivity_estimator=None, check_orthonormality=None, check_tol=None)[source]¶ Bases:
pymor.reductors.basic.StationaryRBReductorReductor for linear
StationaryModelswith affinely decomposed operator and rhs.Note
The reductor
CoerciveRBReductorcan be used for arbitrary coerciveStationaryModelsand offers an improved error estimator with better numerical stability.The only addition is to
StationaryRBReductoris an error estimator, which evaluates the norm of the residual with respect to a given inner product.Parameters
- fom
The
Modelwhich is to be reduced.- RB
VectorArraycontaining the reduced basis on which to project.- product
Inner product for the orthonormalization of
RB, the projection of theOperatorsgiven byvector_ranged_operatorsand for the computation of Riesz representatives of the residual. IfNone, the Euclidean product is used.- coercivity_estimator
Noneor aParameterFunctionalreturning 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.
Methods
assemble_estimator,assemble_estimator_for_subbasisbuild_rom,project_operators,project_operators_to_subbasisextend_basis,reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
h2 module¶
Reductors based on H2-norm.
-
class
pymor.reductors.h2.GenericIRKAReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric IRKA related reductor.
Attributes
-
class
pymor.reductors.h2.IRKAReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.h2.GenericIRKAReductorIterative Rational Krylov Algorithm reductor.
Attributes
-
reduce(rom0_params, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, projection='orth', conv_crit='sigma', compute_errors=False)[source]¶ Reduce using IRKA.
See [GAB08] (Algorithm 4.1) and [ABG10] (Algorithm 1).
Parameters
- rom0_params
Can be:
order of the reduced model (a positive integer),
initial interpolation points (a 1D
NumPy array),dict with
'sigma','b','c'as keys mapping to initial interpolation points (a 1DNumPy array), right tangential directions (VectorArrayfromfom.input_space), and left tangential directions (VectorArrayfromfom.output_space), all of the same length (the order of the reduced model),initial reduced-order model (
LTIModel).
If the order of reduced model is given, initial interpolation data is generated randomly.
- 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'arnoldi': projection matrices are orthogonalized using the Arnoldi process (available only for SISO systems).
- 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
LTIModelmodel.
-
-
class
pymor.reductors.h2.OneSidedIRKAReductor(fom, version, mu=None)[source]¶ Bases:
pymor.reductors.h2.GenericIRKAReductorOne-Sided Iterative Rational Krylov Algorithm reductor.
Parameters
Attributes
-
reduce(rom0_params, 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
- rom0_params
Can be:
order of the reduced model (a positive integer),
initial interpolation points (a 1D
NumPy array),dict with
'sigma','b','c'as keys mapping to initial interpolation points (a 1DNumPy array), right tangential directions (VectorArrayfromfom.input_space), and left tangential directions (VectorArrayfromfom.output_space), all of the same length (the order of the reduced model),initial reduced-order model (
LTIModel).
If the order of reduced model is given, initial interpolation data is generated randomly.
- 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 the current reduced-order model’s poles. Otherwise, only poles in the left half-plane are reflected.- 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
LTIModelmodel.
-
-
class
pymor.reductors.h2.TFIRKAReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.h2.GenericIRKAReductorRealization-independent IRKA reductor.
See [BG12].
Attributes
-
reduce(rom0_params, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, conv_crit='sigma', compute_errors=False)[source]¶ Reduce using TF-IRKA.
Parameters
- rom0_params
Can be:
order of the reduced model (a positive integer),
initial interpolation points (a 1D
NumPy array),dict with
'sigma','b','c'as keys mapping to initial interpolation points (a 1DNumPy array), right tangential directions (VectorArrayfromfom.input_space), and left tangential directions (VectorArrayfromfom.output_space), all of the same length (the order of the reduced model),initial reduced-order model (
LTIModel).
If the order of reduced model is given, initial interpolation data is generated randomly.
- 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 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
- 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
LTIModelmodel.
-
-
class
pymor.reductors.h2.TSIAReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.h2.GenericIRKAReductorTwo-Sided Iteration Algorithm reductor.
Attributes
-
reduce(rom0_params, 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_params
Can be:
order of the reduced model (a positive integer),
initial interpolation points (a 1D
NumPy array),dict with
'sigma','b','c'as keys mapping to initial interpolation points (a 1DNumPy array), right tangential directions (VectorArrayfromfom.input_space), and left tangential directions (VectorArrayfromfom.output_space), all of the same length (the order of the reduced model),initial reduced-order model (
LTIModel).
If the order of reduced model is given, initial interpolation data is generated randomly.
- 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._lti_to_poles_b_c(rom)[source]¶ Compute poles and residues.
Parameters
- rom
Reduced
LTIModel(consisting ofNumpyMatrixOperators).
Returns
- poles
1D
NumPy arrayof poles.- b
VectorArrayfromrom.B.source.- c
VectorArrayfromrom.C.range.
-
pymor.reductors.h2._poles_b_c_to_lti(poles, b, c)[source]¶ Create an
LTIModelfrom poles and residue rank-1 factors.Returns an
LTIModelwith real matrices such that its transfer function is\sum_{i = 1}^r \frac{c_i b_i^T}{s - \lambda_i}
where \lambda_i, b_i, c_i are the poles and residue rank-1 factors.
Parameters
- poles
Sequence of poles.
- b
VectorArrayof right residue rank-1 factors.- c
VectorArrayof left residue rank-1 factors.
Returns
interpolation module¶
-
class
pymor.reductors.interpolation.DelayBHIReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.interpolation.GenericBHIReductorBitangential Hermite interpolation for
LinearDelayModels.Parameters
- fom
The full-order
LinearDelayModelto reduce.- mu
Attributes
-
class
pymor.reductors.interpolation.GenericBHIReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric bitangential Hermite interpolation reductor.
This is a generic reductor for reducing any linear
InputStateOutputModelwith the transfer function which can be written in the generalized coprime factorization H(s) = \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. Interpolation points are assumed to be pairwise distinct.In particular, given interpolation points \sigma_i, right tangential directions b_i, and left tangential directions c_i, for i = 1, 2, \ldots, r, which are closed under conjugation (if \sigma_i is real, then so are b_i and c_i; if \sigma_i is complex, there is \sigma_j such that \sigma_j = \overline{\sigma_i}, b_j = \overline{b_i}, c_j = \overline{c_i}), this reductor finds a transfer function \widehat{H} such that
H(\sigma_i) b_i & = \widehat{H}(\sigma_i) b_i, \\ c_i^T H(\sigma_i) & = c_i^T \widehat{H}(\sigma_i) b_i, \\ c_i^T H'(\sigma_i) b_i & = c_i^T \widehat{H}'(\sigma_i) b_i,
for all i = 1, 2, \ldots, r.
Attributes
-
reduce(sigma, b, c, projection='orth')[source]¶ Bitangential Hermite interpolation.
Parameters
- sigma
Interpolation points (closed under conjugation), sequence of length
r.- b
Right tangential directions,
VectorArrayof lengthrfromself.fom.input_space.- c
Left tangential directions,
VectorArrayof lengthrfromself.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-order model.
-
-
class
pymor.reductors.interpolation.LTIBHIReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.interpolation.GenericBHIReductorBitangential Hermite interpolation for
LTIModels.Attributes
-
reduce(sigma, b, c, projection='orth')[source]¶ Bitangential Hermite interpolation.
Parameters
- sigma
Interpolation points (closed under conjugation), sequence of length
r.- b
Right tangential directions,
VectorArrayof lengthrfromself.fom.input_space.- c
Left tangential directions,
VectorArrayof lengthrfromself.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'arnoldi': projection matrices are orthogonalized using the rational Arnoldi process (available only for SISO systems).
Returns
- rom
Reduced-order model.
-
-
class
pymor.reductors.interpolation.SOBHIReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.interpolation.GenericBHIReductorBitangential Hermite interpolation for
SecondOrderModels.Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
class
pymor.reductors.interpolation.TFBHIReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceLoewner bitangential Hermite interpolation reductor.
See [BG12].
Attributes
-
reduce(sigma, b, c)[source]¶ Realization-independent tangential Hermite interpolation.
Parameters
- sigma
Interpolation points (closed under conjugation), sequence of length
r.- b
Right tangential directions,
VectorArrayfromfom.input_spaceof lengthr.- c
Left tangential directions,
VectorArrayfromfom.output_spaceof lengthr.
Returns
- lti
The reduced-order
LTIModelinterpolating the transfer function offom.
-
parabolic module¶
-
class
pymor.reductors.parabolic.ParabolicRBEstimator(residual, residual_range_dims, initial_residual, initial_residual_range_dims, coercivity_estimator)[source]¶ Bases:
pymor.core.interfaces.ImmutableInterfaceInstantiated by
ParabolicRBReductor.Not to be used directly.
Methods
estimate,restricted_to_subbasisdisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.parabolic.ParabolicRBReductor(fom, RB=None, product=None, coercivity_estimator=None, check_orthonormality=None, check_tol=None)[source]¶ Bases:
pymor.reductors.basic.InstationaryRBReductorReduced Basis Reductor for parabolic equations.
This reductor uses
InstationaryRBReductorfor 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^{-2}(\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
InstationaryModelwhich is to be reduced.- RB
VectorArraycontaining the reduced basis on which to project.- product
The energy inner product
Operatorw.r.t. which the reduction error is estimated andRBis orthonormalized.- coercivity_estimator
Noneor aParameterFunctionalreturning a lower bound C_a(\mu) for the coercivity constant offom.operatorw.r.t.product.
Methods
assemble_estimator,assemble_estimator_for_subbasisbuild_rom,project_operators,project_operators_to_subbasisextend_basis,reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
residual module¶
-
class
pymor.reductors.residual.ImplicitEulerResidualOperator(operator, mass, rhs, dt, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBaseInstantiated by
ImplicitEulerResidualReductor.Methods
apply,projected_to_subbasisapply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,assemble,d_mu,jacobian,pairwise_apply2disable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
H,linear,range,solver_options,source-
apply(U, U_old, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
Parameterfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
-
class
pymor.reductors.residual.ImplicitEulerResidualReductor(RB, operator, mass, dt, rhs=None, product=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceReduced 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_residualusingestimate_image_hierarchical, computes an orthonormal basisresidual_rangeof this range space and then returns the Petrov-Galerkin projectionprojected_riesz_residual == riesz_residual.projected(range_basis=residual_range, source_basis=RB)
of the
riesz_residualoperator. Given reduced basis coefficient vectorsuandu_old, the dual norm of the residual can then be computed asprojected_riesz_residual.apply(u, u_old, mu).l2_norm()
Moreover, a
reconstructmethod is provided such thatresidual_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. IfNone, zero right-hand side is assumed.- RB
VectorArraycontaining a basis of the reduced space onto which to project.- product
Inner product
Operatorw.r.t. which to compute the Riesz representatives.
Methods
reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
-
class
pymor.reductors.residual.NonProjectedImplicitEulerResidualOperator(operator, mass, rhs, dt, product)[source]¶ Bases:
pymor.reductors.residual.ImplicitEulerResidualOperatorInstantiated by
ImplicitEulerResidualReductor.Not to be used directly.
Methods
apply,projected_to_subbasisapply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,assemble,d_mu,jacobian,pairwise_apply2disable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
H,linear,range,solver_options,source-
apply(U, U_old, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
Parameterfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
-
class
pymor.reductors.residual.NonProjectedResidualOperator(operator, rhs, riesz_representatives, product)[source]¶ Bases:
pymor.reductors.residual.ResidualOperatorInstantiated by
ResidualReductor.Not to be used directly.
Methods
apply,projected_to_subbasisapply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,assemble,d_mu,jacobian,pairwise_apply2disable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
H,linear,range,solver_options,source-
apply(U, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
Parameterfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
-
class
pymor.reductors.residual.ResidualOperator(operator, rhs, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBaseInstantiated by
ResidualReductor.Methods
apply,projected_to_subbasisapply2,apply_adjoint,apply_inverse,apply_inverse_adjoint,as_range_array,as_source_array,assemble,d_mu,jacobian,pairwise_apply2disable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
H,linear,range,solver_options,source-
apply(U, mu=None)[source]¶ Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
Parameterfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
-
-
class
pymor.reductors.residual.ResidualReductor(RB, operator, rhs=None, product=None, riesz_representatives=False)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric 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
RBof a subspace of the source space ofoperator, this reductor usesestimate_image_hierarchicalto determine a low-dimensional subspace containing the image of the subspace underresidual(resp.riesz_residual), computes an orthonormal basisresidual_rangefor this range space and then returns the Petrov-Galerkin projectionprojected_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 asprojected_residual.apply(u, mu).l2_norm()
Moreover, a
reconstructmethod is provided such thatresidual_reductor.reconstruct(projected_residual.apply(u, mu)) == residual.apply(RB.lincomb(u), mu)
Parameters
- RB
VectorArraycontaining a basis of the reduced space onto which to project.- operator
See definition of
residual.- rhs
See definition of
residual. IfNone, zero right-hand side is assumed.- product
Inner product
Operatorw.r.t. which to orthonormalize and w.r.t. which to compute the Riesz representatives in caseoperatormaps to functionals.- riesz_representatives
If
Truecompute the Riesz representative of the residual.
Methods
reconstruct,reducedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
sobt module¶
-
class
pymor.reductors.sobt.GenericSOBTpvReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric Second-Order Balanced Truncation position/velocity reductor.
See [RS08].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
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-order
SecondOrderModel.
-
class
pymor.reductors.sobt.SOBTReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceSecond-Order Balanced Truncation reductor.
See [CLVV06].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
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-order
SecondOrderModel.
-
class
pymor.reductors.sobt.SOBTfvReductor(fom, mu=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceFree-velocity Second-Order Balanced Truncation reductor.
See [MS96].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
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-order
SecondOrderModel.
-
class
pymor.reductors.sobt.SOBTpReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.sobt.GenericSOBTpvReductorSecond-Order Balanced Truncation position reductor.
See [RS08].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
class
pymor.reductors.sobt.SOBTpvReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.sobt.GenericSOBTpvReductorSecond-Order Balanced Truncation position-velocity reductor.
See [RS08].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
class
pymor.reductors.sobt.SOBTvReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.sobt.GenericSOBTpvReductorSecond-Order Balanced Truncation velocity reductor.
See [RS08].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
class
pymor.reductors.sobt.SOBTvpReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.sobt.GenericSOBTpvReductorSecond-Order Balanced Truncation velocity-position reductor.
See [RS08].
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
sor_irka module¶
IRKA-type reductor for SecondOrderModels.
-
class
pymor.reductors.sor_irka.SORIRKAReductor(fom, mu=None)[source]¶ Bases:
pymor.reductors.h2.GenericIRKAReductorSOR-IRKA reductor.
Parameters
- fom
The full-order
SecondOrderModelto reduce.- mu
Attributes
-
reduce(rom0_params, tol=0.0001, maxit=100, num_prev=1, force_sigma_in_rhp=False, projection='orth', 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
- rom0_params
Can be:
order of the reduced model (a positive integer),
dict with
'sigma','b','c'as keys mapping to initial interpolation points (a 1DNumPy array), right tangential directions (VectorArrayfromfom.input_space), and left tangential directions (VectorArrayfromfom.output_space), all of the same length (the order of the reduced model),initial reduced-order model (
LTIModel).
If the order of reduced model is given, initial interpolation data is generated randomly.
- 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-order
SecondOrderModel.