pymor.reductors package¶
Submodules¶
basic module¶
-
class
pymor.reductors.basic.GenericRBReductor(d, RB=None, orthogonal_projection=('initial_data', ), product=None)[source]¶ Bases:
pymor.core.interfaces.BasicInterfaceGeneric reduced basis reductor.
Replaces each
Operatorof the givenDiscretizationwith the Galerkin projection onto the span of the given reduced basis.Parameters
- d
- The
Discretizationwhich is to be reduced. - RB
VectorArraycontaining the reduced basis on which to project.- orthogonal_projection
- List of keys in
d.operatorsfor which the correspondingOperatorshould be orthogonally 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
Operatorsgiven byorthogonal_projection.
Methods
GenericRBReductorextend_basis,reconstruct,reduceBasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
BasicInterfacelogger,logging_disabled,name,uid-
extend_basis(U, method='gram_schmidt', pod_modes=1, pod_orthonormalize=True, copy_U=True)[source]¶ Extend basis by new vectors.
Parameters
- U
VectorArraycontaining the new basis vectors.- method
Basis extension method to use. The following methods are available:
trivial: Vectors in Uare appended to the basis. Duplicate vectors in the sense ofalmost_equalare removed.gram_schmidt: New basis vectors are orthonormalized w.r.t. to the old basis using the gram_schmidtalgorithm.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
Trueandmethod == 'pod', re-orthonormalize the new basis vectors obtained by the POD in order to improve numerical accuracy. - copy_U
- If
copy_UisFalse, the new basis vectors might be removed fromU.
Raises
- ExtensionError
- Raised when the selected extension method does not yield a basis of increased dimension.
-
reduce()[source]¶ Perform the reduced basis projection.
Returns
The reduced
Discretization.
-
pymor.reductors.basic.reduce_to_subbasis(d, dim)[source]¶ Further reduce a
Discretizationto the subbasis formed by the firstdimbasis vectors.This is achieved by calling
project_to_subbasisfor each operator of the givenDiscretization.Parameters
- d
- The
Discretizationto further reduce. - dim
- The dimension of the subbasis.
Returns
The further reduced
Discretization.
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
CoerciveRBEstimatorestimate,restricted_to_subbasisImmutableInterfacegenerate_sid,with_,__setattr__BasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementors
-
class
pymor.reductors.coercive.CoerciveRBReductor(d, RB=None, orthogonal_projection=('initial_data', ), product=None, coercivity_estimator=None)[source]¶ Bases:
pymor.reductors.basic.GenericRBReductorReduced Basis reductor for
StationaryDiscretizationswith coercive linear operator.The only addition to
GenericRBReductoris 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].[BEOR14] (1, 2) A. Buhr, C. Engwer, M. Ohlberger, S. Rave, A Numerically Stable A Posteriori Error Estimator for Reduced Basis Approximations of Elliptic Equations, Proceedings of the 11th World Congress on Computational Mechanics, 2014. Parameters
- d
- The
Discretizationwhich is to be reduced. - RB
VectorArraycontaining the reduced basis on which to project.- orthogonal_projection
- List of keys in
d.operatorsfor which the correspondingOperatorshould be orthogonally 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
Operatorsgiven byorthogonal_projectionand 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
CoerciveRBReductorreduceGenericRBReductorextend_basis,reconstructBasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
BasicInterfacelogger,logging_disabled,name,uid
-
class
pymor.reductors.coercive.SimpleCoerciveRBEstimator(estimator_matrix, coercivity_estimator)[source]¶ Bases:
pymor.core.interfaces.ImmutableInterfaceInstantiated by
SimpleCoerciveRBReductor.Not to be used directly.
Methods
SimpleCoerciveRBEstimatorestimate,restricted_to_subbasisImmutableInterfacegenerate_sid,with_,__setattr__BasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementors
-
class
pymor.reductors.coercive.SimpleCoerciveRBReductor(d, RB=None, orthogonal_projection=('initial_data', ), product=None, coercivity_estimator=None)[source]¶ Bases:
pymor.reductors.basic.GenericRBReductorReductor for linear
StationaryDiscretizationswith affinely decomposed operator and rhs.Note
The reductor
CoerciveRBReductorcan be used for arbitrary coerciveStationaryDiscretizationsand offers an improved error estimator with better numerical stability.The only addition is to
GenericRBReductoris an error estimator, which evaluates the norm of the residual with respect to a given inner product.Parameters
- d
- The
Discretizationwhich is to be reduced. - RB
VectorArraycontaining the reduced basis on which to project.- orthogonal_projection
- List of keys in
d.operatorsfor which the correspondingOperatorshould be orthogonally 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
Operatorsgiven byorthogonal_projectionand 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
SimpleCoerciveRBReductorreduceGenericRBReductorextend_basis,reconstructBasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
BasicInterfacelogger,logging_disabled,name,uid
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
ParabolicRBEstimatorestimate,restricted_to_subbasisImmutableInterfacegenerate_sid,with_,__setattr__BasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementors
-
class
pymor.reductors.parabolic.ParabolicRBReductor(d, RB=None, product=None, coercivity_estimator=None)[source]¶ Bases:
pymor.reductors.basic.GenericRBReductorReduced Basis Reductor for parabolic equations.
This reductor uses
GenericRBReductorfor 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}](../_images/math/675f85a40f4b048afec8e48b7a64871dcf8fc730.png)
Here,
denotes the norm induced by the problem’s mass matrix
(e.g. the L^2-norm) and
is an arbitrary energy norm w.r.t.
which the space operator
is coercive, and
is a
lower bound for its coercivity constant. Finally,
denotes
the implicit Euler timestepping residual for the (fixed) time step size
,
where
denotes the mass operator and
the source term.
The dual norm of the residual is computed using the numerically stable projection
from [BEOR14].Warning
The reduced basis
RBis required to be orthonormal w.r.t. the given energy product. If not, the projection of the initial values will be computed incorrectly.[GP05] M. A. Grepl, A. T. Patera, A Posteriori Error Bounds For Reduced-Basis Approximations Of Parametrized Parabolic Partial Differential Equations, M2AN 39(1), 157-181, 2005. [HO08] B. Haasdonk, M. Ohlberger, Reduced basis method for finite volume approximations of parametrized evolution equations, M2AN 42(2), 277-302, 2008. Parameters
- discretization
- The
InstationaryDiscretizationwhich is to be reduced. - RB
VectorArraycontaining the reduced basis on which to project.- product
- The energy inner product
Operatorw.r.t. the reduction error is estimated. RB must be to be orthonomrmal w.r.t. this product! - coercivity_estimator
Noneor aParameterFunctionalreturning a lower bound
for the coercivity constant of discretization.operatorw.r.t.product.
Methods
ParabolicRBReductorreduceGenericRBReductorextend_basis,reconstructBasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
BasicInterfacelogger,logging_disabled,name,uid
residual module¶
-
class
pymor.reductors.residual.ImplicitEulerResidualOperator(operator, mass, functional, dt, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBaseInstantiated by
ImplicitEulerResidualReductor.Methods
-
class
pymor.reductors.residual.ImplicitEulerResidualReductor(RB, operator, mass, dt, functional=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) - functional.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. - functional
- 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
ImplicitEulerResidualReductorreconstruct,reduceBasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
BasicInterfacelogger,logging_disabled,name,uid
-
class
pymor.reductors.residual.NonProjectedImplicitEulerResidualOperator(operator, mass, functional, dt, product)[source]¶ Bases:
pymor.reductors.residual.ImplicitEulerResidualOperatorInstantiated by
ImplicitEulerResidualReductor.Not to be used directly.
Methods
-
class
pymor.reductors.residual.NonProjectedResidualOperator(operator, rhs, rhs_is_functional, product)[source]¶ Bases:
pymor.reductors.residual.ResidualOperatorInstantiated by
ResidualReductor.Not to be used directly.
Methods
-
class
pymor.reductors.residual.ResidualOperator(operator, rhs, rhs_is_functional=True, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBaseInstantiated by
ResidualReductor.Methods
-
class
pymor.reductors.residual.ResidualReductor(RB, operator, rhs=None, product=None)[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_vector(mu)
When the rhs is a functional we are interested in the Riesz representative of the residual:
residual.apply(U, mu) == product.apply_inverse(operator.apply(U, mu) - rhs.as_vector(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 caserhsis a functional.
Methods
ResidualReductorreconstruct,reduceBasicInterfacedisable_logging,enable_logging,has_interface_name,implementor_names,implementorsAttributes
BasicInterfacelogger,logging_disabled,name,uid