pymor.operators package¶
Submodules¶
basic module¶
-
class
pymor.operators.basic.
OperatorBase
[source]¶ Bases:
pymor.operators.interfaces.OperatorInterface
Base class for
Operators
providing some default implementations.When implementing a new operator, it is usually advisable to derive from this class.
Methods
-
class
pymor.operators.basic.
ProjectedOperator
(operator, range_basis, source_basis, product=None, solver_options=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Generic
Operator
representing the projection of anOperator
to a subspace.This operator is implemented as the concatenation of the linear combination with
source_basis
, application of the originaloperator
and projection ontorange_basis
. As such, this operator can be used to obtain a reduced basis projection of any givenOperator
. However, no offline/online decomposition is performed, so this operator is mainly useful for testing before implementing offline/online decomposition for a specific application.This operator is instantiated in
pymor.algorithms.projection.project
as a default implementation for parametric or nonlinear operators.Parameters
- operator
- The
Operator
to project. - source_basis
- See
pymor.algorithms.projection.project
. - range_basis
- See
pymor.algorithms.projection.project
. - product
- See
pymor.algorithms.projection.project
.
Methods
block module¶
-
class
pymor.operators.block.
BlockDiagonalOperator
(blocks)[source]¶ Bases:
pymor.operators.block.BlockOperator
Block diagonal
Operator
of arbitraryOperators
.This is a specialization of
BlockOperator
for the block diagonal case.Methods
-
class
pymor.operators.block.
BlockOperator
(blocks)[source]¶ Bases:
pymor.operators.basic.OperatorBase
A matrix of arbitrary
Operators
.This operator can be
applied
to a compatibleBlockVectorArrays
.Parameters
- blocks
- Two-dimensional array-like where each entry is an
Operator
orNone
.
Methods
Attributes
cg module¶
This module provides some operators for continuous finite element discretizations.
-
class
pymor.operators.cg.
AdvectionOperatorP1
(grid, boundary_info, advection_function=None, advection_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear advection
Operator
for linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
v
is vector-valued. The current implementation works in one and two dimensions, but can be trivially extended to arbitrary dimensions.Parameters
- grid
- The
Grid
for which to assemble the operator. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- advection_function
- The
Function
v(x)
withshape_range = (grid.dim, )
. IfNone
, constant one is assumed. - advection_constant
- The constant
c
. IfNone
,c
is set to one. - dirichlet_clear_columns
- If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero. - dirichlet_clear_diag
- If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one. - name
- Name of the operator.
Methods
-
class
pymor.operators.cg.
AdvectionOperatorQ1
(grid, boundary_info, advection_function=None, advection_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear advection
Operator
for bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
v
has to be vector-valued. The current implementation works in two dimensions, but can be trivially extended to arbitrary dimensions.Parameters
- grid
- The
Grid
for which to assemble the operator. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- advection_function
- The
Function
v(x)
withshape_range = (grid.dim, )
. IfNone
, constant one is assumed. - advection_constant
- The constant
c
. IfNone
,c
is set to one. - dirichlet_clear_columns
- If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero. - dirichlet_clear_diag
- If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one. - name
- Name of the operator.
Methods
-
class
pymor.operators.cg.
DiffusionOperatorP1
(grid, boundary_info, diffusion_function=None, diffusion_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Diffusion
Operator
for linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
d
can be scalar- or matrix-valued. The current implementation works in one and two dimensions, but can be trivially extended to arbitrary dimensions.Parameters
- grid
- The
Grid
for which to assemble the operator. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- diffusion_function
- The
Function
d(x)
withshape_range == ()
orshape_range = (grid.dim, grid.dim)
. IfNone
, constant one is assumed. - diffusion_constant
- The constant
c
. IfNone
,c
is set to one. - dirichlet_clear_columns
- If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero. - dirichlet_clear_diag
- If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one. - name
- Name of the operator.
Methods
-
class
pymor.operators.cg.
DiffusionOperatorQ1
(grid, boundary_info, diffusion_function=None, diffusion_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Diffusion
Operator
for bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
d
can be scalar- or matrix-valued. The current implementation works in two dimensions, but can be trivially extended to arbitrary dimensions.Parameters
- grid
- The
Grid
for which to assemble the operator. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- diffusion_function
- The
Function
d(x)
withshape_range == ()
orshape_range = (grid.dim, grid.dim)
. IfNone
, constant one is assumed. - diffusion_constant
- The constant
c
. IfNone
,c
is set to one. - dirichlet_clear_columns
- If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero. - dirichlet_clear_diag
- If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one. - name
- Name of the operator.
Methods
-
class
pymor.operators.cg.
InterpolationOperator
(grid, function)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Vector-like Lagrange interpolation
Operator
for continuous finite element spaces.Methods
-
class
pymor.operators.cg.
L2ProductFunctionalP1
(grid, function, boundary_info=None, dirichlet_data=None, neumann_data=None, robin_data=None, order=2, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear finite element
Functional
representing the inner product with an L2-Function
.Boundary treatment can be performed by providing
boundary_info
anddirichlet_data
, in which case the DOFs corresponding to Dirichlet boundaries are set to the values provided bydirichlet_data
. Neumann boundaries are handled by providing aneumann_data
function, Robin boundaries by providing arobin_data
tuple.The current implementation works in one and two dimensions, but can be trivially extended to arbitrary dimensions.
Parameters
- grid
Grid
for which to assemble the functional.- function
- The
Function
with which to take the inner product. - boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries orNone
. IfNone
, no boundary treatment is performed.- dirichlet_data
Function
providing the Dirichlet boundary values. IfNone
, constant-zero boundary is assumed.- neumann_data
Function
providing the Neumann boundary values. IfNone
, constant-zero is assumed.- robin_data
- Tuple of two
Functions
providing the Robin parameter and boundary values, seeRobinBoundaryOperator
. IfNone
, constant-zero for both functions is assumed. - order
- Order of the Gauss quadrature to use for numerical integration.
- name
- The name of the functional.
Methods
-
class
pymor.operators.cg.
L2ProductFunctionalQ1
(grid, function, boundary_info=None, dirichlet_data=None, neumann_data=None, robin_data=None, order=2, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Bilinear finite element
Functional
representing the inner product with an L2-Function
.Boundary treatment can be performed by providing
boundary_info
anddirichlet_data
, in which case the DOFs corresponding to Dirichlet boundaries are set to the values provided bydirichlet_data
. Neumann boundaries are handled by providing aneumann_data
function, Robin boundaries by providing arobin_data
tuple.The current implementation works in two dimensions, but can be trivially extended to arbitrary dimensions.
Parameters
- grid
Grid
for which to assemble the functional.- function
- The
Function
with which to take the inner product. - boundary_info
BoundaryInfo
determining the Dirichlet boundaries orNone
. IfNone
, no boundary treatment is performed.- dirichlet_data
Function
providing the Dirichlet boundary values. IfNone
, constant-zero boundary is assumed.- neumann_data
Function
providing the Neumann boundary values. IfNone
, constant-zero is assumed.- robin_data
- Tuple of two
Functions
providing the Robin parameter and boundary values, seeRobinBoundaryOperator
. IfNone
, constant-zero for both functions is assumed. - order
- Order of the Gauss quadrature to use for numerical integration.
- name
- The name of the functional.
Methods
-
class
pymor.operators.cg.
L2ProductP1
(grid, boundary_info, dirichlet_clear_rows=True, dirichlet_clear_columns=False, dirichlet_clear_diag=False, coefficient_function=None, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Operator
representing the L2-product between linear finite element functions.The current implementation works in one and two dimensions, but can be trivially extended to arbitrary dimensions.
Parameters
- grid
- The
Grid
for which to assemble the product. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- dirichlet_clear_rows
- If
True
, set the rows of the system matrix corresponding to Dirichlet boundary DOFs to zero. - dirichlet_clear_columns
- If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero. - dirichlet_clear_diag
- If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise, if eitherdirichlet_clear_rows
ordirichlet_clear_columns
isTrue
, the diagonal entries are set to one. - coefficient_function
- Coefficient
Function
for product withshape_range == ()
. IfNone
, constant one is assumed. - name
- The name of the product.
Methods
-
class
pymor.operators.cg.
L2ProductQ1
(grid, boundary_info, dirichlet_clear_rows=True, dirichlet_clear_columns=False, dirichlet_clear_diag=False, coefficient_function=None, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Operator
representing the L2-product between bilinear finite element functions.The current implementation works in two dimensions, but can be trivially extended to arbitrary dimensions.
Parameters
- grid
- The
Grid
for which to assemble the product. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- dirichlet_clear_rows
- If
True
, set the rows of the system matrix corresponding to Dirichlet boundary DOFs to zero. - dirichlet_clear_columns
- If
True
, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero. - dirichlet_clear_diag
- If
True
, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise, if eitherdirichlet_clear_rows
ordirichlet_clear_columns
isTrue
, the diagonal entries are set to one. - coefficient_function
- Coefficient
Function
for product withshape_range == ()
. IfNone
, constant one is assumed. - name
- The name of the product.
Methods
-
class
pymor.operators.cg.
RobinBoundaryOperator
(grid, boundary_info, robin_data=None, order=2, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Robin boundary
Operator
for linear finite elements.The operator represents the contribution of Robin boundary conditions to the stiffness matrix, where the boundary condition is supposed to be given in the form
-[ d(x) ∇u(x) ] ⋅ n(x) = c(x) (u(x) - g(x))
d
andn
are the diffusion function (seeDiffusionOperatorP1
) and the unit outer normal inx
, whilec
is the (scalar) Robin parameter function andg
is the (also scalar) Robin boundary value function.Parameters
- grid
- The
Grid
over which to assemble the operator. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- robin_data
- Tuple providing two
Functions
that represent the Robin parameter and boundary value function. IfNone
, the resulting operator is zero. - name
- Name of the operator.
Methods
constructions module¶
Module containing some constructions to obtain new operators from old ones.
-
class
pymor.operators.constructions.
AdjointOperator
(operator, source_product=None, range_product=None, name=None, with_apply_inverse=True, solver_options=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Represents the adjoint of a given linear
Operator
.For a linear
Operator
op
the adjointop^*
ofop
is given by:(op^*(v), u)_s = (v, op(u))_r,
where
( , )_s
and( , )_r
denote the inner products on the source and range space ofop
. If two products are given byP_s
andP_r
, then:op^*(v) = P_s^(-1) o op.T o P_r,
Thus, if
( , )_s
and( , )_r
are the Euclidean inner products,op^*v
is simply given by applycation of the :attr:transpose <pymor.operators.interface.OperatorInterface.T>`Operator
.Parameters
- operator
- The
Operator
of which the adjoint is formed. - source_product
- If not
None
, inner productOperator
for the sourceVectorSpaceInterface
w.r.t. which to take the adjoint. - range_product
- If not
None
, inner productOperator
for the rangeVectorSpaceInterface
w.r.t. which to take the adjoint. - name
- If not
None
, name of the operator. - with_apply_inverse
- If
True
, provide ownapply_inverse
andapply_inverse_transpose
implementations by calling these methods on the givenoperator
. (Is set toFalse
in the default implementation of andapply_inverse_transpose
.) - solver_options
- When
with_apply_inverse
isFalse
, thesolver_options
to use for theapply_inverse
default implementation.
Methods
-
class
pymor.operators.constructions.
AffineOperator
(operator, name=None)[source]¶ Bases:
pymor.operators.constructions.ProxyOperator
Decompose an affine
Operator
into affine_shift and linear_part.Methods
-
class
pymor.operators.constructions.
ComponentProjection
(components, source, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Operator
representing the projection of aVectorArray
on some of its components.Parameters
- components
- List or 1D
NumPy array
of the indices of the vectorcomponents
that ar to be extracted by the operator. - source
- Source
VectorSpaceInterface
of the operator. - name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
Concatenation
(second, first, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Operator
representing the concatenation of twoOperators
.Parameters
Methods
-
class
pymor.operators.constructions.
ConstantOperator
(value, source, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
A constant
Operator
always returning the same vector.Parameters
- value
- A
VectorArray
of length 1 containing the vector which is returned by the operator. - source
- Source
VectorSpaceInterface
of the operator. - name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
FixedParameterOperator
(operator, mu=None, name=None)[source]¶ Bases:
pymor.operators.constructions.ProxyOperator
Makes an
Operator
Parameter
-independent by setting a fixedParameter
.Parameters
Methods
-
class
pymor.operators.constructions.
IdentityOperator
(space, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
The identity
Operator
.In other words:
op.apply(U) == U
Parameters
- space
- The
VectorSpaceInterface
the operator acts on. - name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
InducedNorm
(product, raise_negative, tol, name)[source]¶ Bases:
pymor.core.interfaces.ImmutableInterface
,pymor.parameters.base.Parametric
Instantiated by
induced_norm
. Do not use directly.
-
class
pymor.operators.constructions.
InverseOperator
(operator, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Represents the inverse of a given
Operator
.Parameters
- operator
- The
Operator
of which the inverse is formed. - name
- If not
None
, name of the operator.
Methods
-
class
pymor.operators.constructions.
InverseTransposeOperator
(operator, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Represents the inverse transpose of a given
Operator
.Parameters
- operator
- The
Operator
of which the inverse transpose is formed. - name
- If not
None
, name of the operator.
Methods
-
class
pymor.operators.constructions.
LincombOperator
(operators, coefficients, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Linear combination of arbitrary
Operators
.This
Operator
represents a (possiblyParameter
dependent) linear combination of a given list ofOperators
.Parameters
- operators
- List of
Operators
whose linear combination is formed. - coefficients
- A list of linear coefficients. A linear coefficient can
either be a fixed number or a
ParameterFunctional
. - name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
LinearOperator
(operator, name=None)[source]¶ Bases:
pymor.operators.constructions.ProxyOperator
Mark the wrapped
Operator
to be linear.Methods
-
class
pymor.operators.constructions.
ProxyOperator
(operator, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Forwards all interface calls to given
Operator
.Mainly useful as base class for other
Operator
implementations.Parameters
- operator
- The
Operator
to wrap. - name
- Name of the wrapping operator.
Methods
-
class
pymor.operators.constructions.
SelectionOperator
(operators, parameter_functional, boundaries, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
An
Operator
selected from a list ofOperators
.operators[i]
is used ifparameter_functional(mu)
is less or equal thanboundaries[i]
and greater thanboundaries[i-1]
:-infty ------- boundaries[i] ---------- boundaries[i+1] ------- infty | | --- operators[i] ---|---- operators[i+1] ----|---- operators[i+2] | |
Parameters
- operators
- List of
Operators
from which oneOperator
is selected based on the givenParameter
. - parameter_functional
- The
ParameterFunctional
used for the selection of oneOperator
. - boundaries
- The interval boundaries as defined above.
- name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
VectorArrayOperator
(array, transposed=False, space_id=None, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Wraps a
VectorArray
as anOperator
.If
transposed
isFalse
, the operator maps fromNumpyVectorSpace(len(array))
toarray.space
by forming linear combinations of the vectors in the array with given coefficient arrays.If
transposed == True
, the operator maps fromarray.space
toNumpyVectorSpace(len(array))
by forming the inner products of the argument with the vectors in the given array.Parameters
- array
- The
VectorArray
which is to be treated as an operator. - transposed
- See description above.
- name
- The name of the operator.
Methods
-
class
pymor.operators.constructions.
VectorFunctional
(vector, product=None, name=None)[source]¶ Bases:
pymor.operators.constructions.VectorArrayOperator
Wrap a vector as a linear
Functional
.Given a vector
v
of dimensiond
, this class represents the functionalf: R^d ----> R^1 u |---> (u, v)
where
( , )
denotes the inner product given byproduct
.In particular, if
product
isNone
VectorFunctional(vector).as_source_array() == vector.
If
product
is not none, we obtainVectorFunctional(vector).as_source_array() == product.apply(vector).
Parameters
- vector
VectorArray
of length 1 containing the vectorv
.- product
Operator
representing the scalar product to use.- name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
VectorOperator
(vector, name=None)[source]¶ Bases:
pymor.operators.constructions.VectorArrayOperator
Wrap a vector as a vector-like
Operator
.Given a vector
v
of dimensiond
, this class represents the operatorop: R^1 ----> R^d x |---> x⋅v
In particular:
VectorOperator(vector).as_range_array() == vector
Parameters
- vector
VectorArray
of length 1 containing the vectorv
.- name
- Name of the operator.
Methods
-
class
pymor.operators.constructions.
ZeroOperator
(source, range, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
The
Operator
which maps every vector to zero.Parameters
- source
- Source
VectorSpaceInterface
of the operator. - range
- Range
VectorSpaceInterface
of the operator. - name
- Name of the operator.
Methods
-
pymor.operators.constructions.
induced_norm
(product, raise_negative=True, tol=1e-10, name=None)[source]¶ Obtain induced norm of an inner product.
The norm of the vectors in a
VectorArray
U is calculated by callingproduct.pairwise_apply2(U, U, mu=mu).
In addition, negative norm squares of absolute value smaller than
tol
are clipped to0
. Ifraise_negative
isTrue
, aValueError
exception is raised if there are negative norm squares of absolute value larger thantol
.Parameters
- product
- The inner product
Operator
for which the norm is to be calculated. - raise_negative
- If
True
, raise an exception if calculated norm is negative. - tol
- See above.
- name
- optional, if None product’s name is used
Returns
- norm
- A function
norm(U, mu=None)
taking aVectorArray
U
as input together with theParameter
mu
which is passed to the product.
Defaults
raise_negative, tol (see
pymor.core.defaults
)
ei module¶
-
class
pymor.operators.ei.
EmpiricalInterpolatedOperator
(operator, interpolation_dofs, collateral_basis, triangular, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Interpolate an
Operator
using Empirical Operator Interpolation.Let
L
be anOperator
,0 <= c_1, ..., c_M < L.range.dim
indices of interpolation DOFs and letb_1, ..., b_M in R^(L.range.dim)
be collateral basis vectors. If moreoverψ_j(U)
denotes the j-th component ofU
, the empirical interpolationL_EI
ofL
w.r.t. the given data is given by| M | L_EI(U, μ) = ∑ b_i⋅λ_i such that | i=1 | | ψ_(c_i)(L_EI(U, μ)) = ψ_(c_i)(L(U, μ)) for i=0,...,M
Since the original operator only has to be evaluated at the given interpolation DOFs,
EmpiricalInterpolatedOperator
callsrestricted
to obtain a restricted version of the operator which is used to quickly obtain the required evaluations. If therestricted
method, is not implemented, the full operator will be evaluated (which will lead to the same result, but without any speedup).The interpolation DOFs and the collateral basis can be generated using the algorithms provided in the
pymor.algorithms.ei
module.Parameters
- operator
- The
Operator
to interpolate. - interpolation_dofs
- List or 1D
NumPy array
of the interpolation DOFsc_1, ..., c_M
. - collateral_basis
VectorArray
containing the collateral basisb_1, ..., b_M
.- triangular
- If
True
, assume that ψ_(c_i)(b_j) = 0 for i < j, which means that the interpolation matrix is triangular. - name
- Name of the operator.
Methods
-
class
pymor.operators.ei.
ProjectedEmpiciralInterpolatedOperator
(restricted_operator, interpolation_matrix, source_basis_dofs, projected_collateral_basis, triangular, source_id, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
A projected
EmpiricalInterpolatedOperator
.Methods
fv module¶
This module provides some operators for finite volume discretizations.
-
class
pymor.operators.fv.
DiffusionOperator
(grid, boundary_info, diffusion_function=None, diffusion_constant=None, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Finite Volume Diffusion
Operator
.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
Parameters
- grid
- The
Grid
over which to assemble the operator. - boundary_info
BoundaryInfo
for the treatment of Dirichlet boundary conditions.- diffusion_function
- The scalar-valued
Function
d(x)
. IfNone
, constant one is assumed. - diffusion_constant
- The constant
c
. IfNone
,c
is set to one. - name
- Name of the operator.
Methods
-
class
pymor.operators.fv.
EngquistOsherFlux
(flux, flux_derivative, gausspoints=5, intervals=1)[source]¶ Bases:
pymor.operators.fv.NumericalConvectiveFluxInterface
Engquist-Osher numerical flux.
If
f
is the analytical flux, andf'
its derivative, the Engquist-Osher flux is given by:F(U_in, U_out, normal, vol) = vol * [c^+(U_in, normal) + c^-(U_out, normal)] U_in c^+(U_in, normal) = f(0)⋅normal + ∫ max(f'(s)⋅normal, 0) ds s=0 U_out c^-(U_out, normal) = ∫ min(f'(s)⋅normal, 0) ds s=0
Parameters
Methods
EngquistOsherFlux
evaluate_stage1
,evaluate_stage2
ImmutableInterface
generate_sid
,with_
,__setattr__
BasicInterface
disable_logging
,enable_logging
,has_interface_name
,implementor_names
,implementors
Parametric
build_parameter_type
,parse_parameter
,strip_parameter
-
class
pymor.operators.fv.
L2Product
(grid, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Operator
representing the L2-product between finite volume functions.Parameters
- grid
- The
Grid
for which to assemble the product. - name
- The name of the product.
Methods
-
class
pymor.operators.fv.
L2ProductFunctional
(grid, function=None, boundary_info=None, dirichlet_data=None, diffusion_function=None, diffusion_constant=None, neumann_data=None, order=1, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Finite volume
Functional
representing the inner product with an L2-Function
.Additionally, boundary conditions can be enforced by providing
dirichlet_data
andneumann_data
functions.Parameters
- grid
Grid
for which to assemble the functional.- function
- The
Function
with which to take the inner product orNone
. - boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries orNone
. IfNone
, no boundary treatment is performed.- dirichlet_data
Function
providing the Dirichlet boundary values. IfNone
, constant-zero boundary is assumed.- diffusion_function
- See
DiffusionOperator
. Has to be specified in casedirichlet_data
is given. - diffusion_constant
- See
DiffusionOperator
. Has to be specified in casedirichlet_data
is given. - neumann_data
Function
providing the Neumann boundary values. IfNone
, constant-zero is assumed.- order
- Order of the Gauss quadrature to use for numerical integration.
- name
- The name of the functional.
Methods
-
class
pymor.operators.fv.
LaxFriedrichsFlux
(flux, lxf_lambda=1.0)[source]¶ Bases:
pymor.operators.fv.NumericalConvectiveFluxInterface
Lax-Friedrichs numerical flux.
If
f
is the analytical flux, the Lax-Friedrichs fluxF
is given by:F(U_in, U_out, normal, vol) = vol * [normal⋅(f(U_in) + f(U_out))/2 + (U_in - U_out)/(2*λ)]
Parameters
- flux
Function
defining the analytical fluxf
.- lxf_lambda
- The stabilization parameter
λ
.
Methods
LaxFriedrichsFlux
evaluate_stage1
,evaluate_stage2
ImmutableInterface
generate_sid
,with_
,__setattr__
BasicInterface
disable_logging
,enable_logging
,has_interface_name
,implementor_names
,implementors
Parametric
build_parameter_type
,parse_parameter
,strip_parameter
-
class
pymor.operators.fv.
LinearAdvectionLaxFriedrichs
(grid, boundary_info, velocity_field, lxf_lambda=1.0, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Linear advection finite Volume
Operator
using Lax-Friedrichs flux.The operator is of the form
L(u, mu)(x) = ∇ ⋅ (v(x, mu)⋅u(x))
See
LaxFriedrichsFlux
for the definition of the Lax-Friedrichs flux.Parameters
- grid
Grid
over which to assemble the operator.- boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries.- velocity_field
Function
defining the velocity fieldv
.- lxf_lambda
- The stabilization parameter
λ
. - name
- The name of the operator.
Methods
-
class
pymor.operators.fv.
NonlinearAdvectionOperator
(grid, boundary_info, numerical_flux, dirichlet_data=None, solver_options=None, space_id='STATE', name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Nonlinear finite volume advection
Operator
.The operator is of the form
L(u, mu)(x) = ∇ ⋅ f(u(x), mu)
Parameters
- grid
Grid
for which to evaluate the operator.- boundary_info
BoundaryInfo
determining the Dirichlet and Neumann boundaries.- numerical_flux
- The
NumericalConvectiveFlux
to use. - dirichlet_data
Function
providing the Dirichlet boundary values. IfNone
, constant-zero boundary is assumed.- name
- The name of the operator.
Methods
-
class
pymor.operators.fv.
NonlinearReactionOperator
(grid, reaction_function, reaction_function_derivative=None, space_id='STATE', name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Methods
-
class
pymor.operators.fv.
NumericalConvectiveFluxInterface
[source]¶ Bases:
pymor.core.interfaces.ImmutableInterface
,pymor.parameters.base.Parametric
Interface for numerical convective fluxes for finite volume schemes.
Numerical fluxes defined by this interfaces are functions of the form
F(U_inner, U_outer, unit_outer_normal, edge_volume, mu)
.- The flux evaluation is vectorized and happens in two stages:
evaluate_stage1
receives aNumPy array
U
of all values which appear asU_inner
orU_outer
for all edges the flux shall be evaluated at and returns atuple
ofNumPy arrays
each of the same length asU
.evaluate_stage2
receives the reorderedstage1_data
for each edge as well as the unit outer normal and the volume of the edges.stage1_data
is given as follows: IfR_l
isl
-th entry of thetuple
returned byevaluate_stage1
, thel
-th entryD_l
of of thestage1_data
tuple has the shape(num_edges, 2) + R_l.shape[1:]
. If for edgek
the valuesU_inner
andU_outer
are thei
-th andj
-th value in theU
array provided toevaluate_stage1
, we haveD_l[k, 0] == R_l[i], D_l[k, 1] == R_l[j].
evaluate_stage2
returns aNumPy array
of the flux evaluations for each edge.
Methods
-
class
pymor.operators.fv.
ReactionOperator
(grid, reaction_coefficient, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Finite Volume reaction
Operator
.The operator is of the form
L(u, mu)(x) = c(x, mu)⋅u(x)
Parameters
- grid
- The
Grid
for which to assemble the operator. - reaction_coefficient
- The function ‘c’
- name
- The name of the operator.
Methods
-
class
pymor.operators.fv.
SimplifiedEngquistOsherFlux
(flux, flux_derivative)[source]¶ Bases:
pymor.operators.fv.NumericalConvectiveFluxInterface
Engquist-Osher numerical flux. Simplified Implementation for special case.
For the definition of the Engquist-Osher flux see
EngquistOsherFlux
. This class provides a faster and more accurate implementation for the special case thatf(0) == 0
and the derivative off
only changes sign at0
.Parameters
Methods
-
pymor.operators.fv.
nonlinear_advection_engquist_osher_operator
(grid, boundary_info, flux, flux_derivative, gausspoints=5, intervals=1, dirichlet_data=None, solver_options=None, name=None)[source]¶ Instantiate a
NonlinearAdvectionOperator
usingEngquistOsherFlux
.
-
pymor.operators.fv.
nonlinear_advection_lax_friedrichs_operator
(grid, boundary_info, flux, lxf_lambda=1.0, dirichlet_data=None, solver_options=None, name=None)[source]¶ Instantiate a
NonlinearAdvectionOperator
usingLaxFriedrichsFlux
.
-
pymor.operators.fv.
nonlinear_advection_simplified_engquist_osher_operator
(grid, boundary_info, flux, flux_derivative, dirichlet_data=None, solver_options=None, name=None)[source]¶ Instantiate a
NonlinearAdvectionOperator
usingSimplifiedEngquistOsherFlux
.
interfaces module¶
-
class
pymor.operators.interfaces.
OperatorInterface
[source]¶ Bases:
pymor.core.interfaces.ImmutableInterface
,pymor.parameters.base.Parametric
Interface for
Parameter
dependent discrete operators.An operator in pyMOR is simply a mapping which for any given
Parameter
maps vectors from its sourceVectorSpaceInterface
to vectors in its rangeVectorSpaceInterface
.Note that there is no special distinction between functionals and operators in pyMOR. A functional is simply an operator with
NumpyVectorSpace
(1)
as its rangeVectorSpaceInterface
.Methods
Attributes
-
solver_options
¶ If not
None
, a dict which can contain the following keys:‘inverse’: solver options used for apply_inverse
‘inverse_transpose’: solver options used for apply_inverse_transpose
‘jacobian’: solver options for the operators returned by jacobian
(has no effect for linear operators)If
solver_options
isNone
or a dict entry is missing orNone
, default options are used. The interpretation of the given solver options is up to the operator at hand. In general, values insolver_options
should either be strings (indicating a solver type) or dicts of options, usually with an entry'type'
which specifies the solver type to use and further items which configure this solver.
-
linear
¶ True
if the operator is linear.
-
source
¶ The source
VectorSpaceInterface
.
-
range
¶ The range
VectorSpaceInterface
.
-
T
¶ The transpose operator, i.e.
self.T.apply(V, mu) == self.apply_transpose(V, mu)
for all V, mu.
-
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.
-
apply2
(V, U, mu=None)[source]¶ Treat the operator as a 2-form by computing
V.dot(self.apply(U))
.If the operator is a linear operator given by multiplication with a matrix M, then
apply2
is given as:op.apply2(V, U) = V^T*M*U.
Parameters
- V
VectorArray
of the left arguments V.- U
VectorArray
of the right right arguments U.- mu
- The
Parameter
for which to evaluate the operator.
Returns
A
NumPy array
with shape(len(V), len(U))
containing the 2-form evaluations.
-
apply_inverse
(V, mu=None, least_squares=False)[source]¶ Apply the inverse operator.
Parameters
- V
VectorArray
of vectors to which the inverse operator is applied.- mu
- The
Parameter
for which to evaluate the inverse operator. - least_squares
If
True
, solve the least squares problem:u = argmin ||op(u) - v||_2.
Since for an invertible operator the least squares solution agrees with the result of the application of the inverse operator, setting this option should, in general, have no effect on the result for those operators. However, note that when no appropriate
solver_options
are set for the operator, most implementations will choose a least squares solver by default which may be undesirable.
Returns
VectorArray
of the inverse operator evaluations.Raises
- InversionError
- The operator could not be inverted.
-
apply_inverse_transpose
(U, mu=None, least_squares=False)[source]¶ Apply the inverse transpose operator.
Parameters
- U
VectorArray
of vectors to which the inverse transpose operator is applied.- mu
- The
Parameter
for which to evaluate the inverse transpose operator. - least_squares
If
True
, solve the least squares problem:v = argmin ||op*(v) - u||_2.
Since for an invertible operator the least squares solution agrees with the result of the application of the inverse operator, setting this option should, in general, have no effect on the result for those operators. However, note that when no appropriate
solver_options
are set for the operator, most operator implementations will choose a least squares solver by default which may be undesirable.
Returns
VectorArray
of the inverse transpose operator evaluations.Raises
- InversionError
- The operator could not be inverted.
-
apply_transpose
(V, mu=None)[source]¶ Apply the transpose operator.
For any given linear
Operator
op
,Parameter
mu
andVectorArrays
U
,V
in thesource
resp.range
we have:op.apply_transpose(V, mu).dot(U) == V.dot(op.apply(U, mu))
Thus, when
op
is represented by a matrixM
,apply_transpose
is given by left-multplication ofM
withV
.Parameters
- V
VectorArray
of vectors to which the transpose operator is applied.- mu
- The
Parameter
for which to apply the transpose operator.
Returns
VectorArray
of the transpose operator evaluations.
-
as_range_array
(mu=None)[source]¶ Return a
VectorArray
representation of the operator in its range space.In the case of a linear operator with
NumpyVectorSpace
assource
, this method returns for everyParameter
mu
aVectorArray
V
in the operator’srange
, such thatV.lincomb(U.data) == self.apply(U, mu)
for all
VectorArrays
U
.Parameters
- mu
- The
Parameter
for which to return theVectorArray
representation.
Returns
- V
- The
VectorArray
defined above.
-
as_source_array
(mu=None)[source]¶ Return a
VectorArray
representation of the operator in its source space.In the case of a linear operator with
NumpyVectorSpace
asrange
, this method returns for everyParameter
mu
aVectorArray
V
in the operator’ssource
, such thatself.range.make_array(V.dot(U)) == self.apply(U, mu)
for all
VectorArrays
U
.Parameters
- mu
- The
Parameter
for which to return theVectorArray
representation.
Returns
- V
- The
VectorArray
defined above.
-
as_vector
(mu=None, *, space=None)[source]¶ Return a vector representation of a linear functional or vector operator.
Depending on the operator’s
source
andrange
, this method is equivalent to callingas_range_array
oras_source_array
respectively. The resultingVectorArray
is required to have length 1.Note that in case both
source
andrange
are one-dimensionalNumpyVectorSpaces
but with differentids
, it is impossible to determine which space to choose. In this case, the desired space has to be specified via thespace
parameter.Parameters
- mu
- The
Parameter
for which to return the vector representation. - space
- See above.
Returns
- V
VectorArray
of length 1 containing the vector representation.
-
assemble
(mu=None)[source]¶ Assemble the operator for a given parameter.
The result of the method strongly depends on the given operator. For instance, a matrix-based operator will assemble its matrix, a
LincombOperator
will try to form the linear combination of its operators, whereas an arbitrary operator might simply return aFixedParameterOperator
. The only assured property of the assembled operator is that it no longer depends on aParameter
.Parameters
- mu
- The
Parameter
for which to assemble the operator.
Returns
Parameter-independent, assembled
Operator
.
-
assemble_lincomb
(operators, coefficients, solver_options=None, name=None)[source]¶ Try to assemble a linear combination of the given operators.
This method is called in the
assemble
method ofLincombOperator
on the first of its operator. If an assembly of the given linear combination is possible, e.g. the linear combination of the system matrices of the operators can be formed, then the assembled operator is returned. Otherwise, the method returnsNone
to indicate that assembly is not possible.Parameters
- operators
- List of
Operators
whose linear combination is formed. - coefficients
- List of the corresponding linear coefficients.
- solver_options
solver_options
for the assembled operator.- name
- Name of the assembled operator.
Returns
The assembled
Operator
if assembly is possible, otherwiseNone
.
-
jacobian
(U, mu=None)[source]¶ Return the operator’s Jacobian as a new
Operator
.Parameters
- U
- Length 1
VectorArray
containing the vector for which to compute the Jacobian. - mu
- The
Parameter
for which to compute the Jacobian.
Returns
Linear
Operator
representing the Jacobian.
-
pairwise_apply2
(V, U, mu=None)[source]¶ Treat the operator as a 2-form by computing
V.dot(self.apply(U))
.Same as
OperatorInterface.apply2
, except that vectors fromV
andU
are applied in pairs.Parameters
- V
VectorArray
of the left arguments V.- U
VectorArray
of the right right arguments U.- mu
- The
Parameter
for which to evaluate the operator.
Returns
A
NumPy array
with shape(len(V),) == (len(U),)
containing the 2-form evaluations.
-
restricted
(dofs)[source]¶ Restrict the operator range to a given set of degrees of freedom.
This method returns a restricted version
restricted_op
of the operator along with an arraysource_dofs
such that for anyVectorArray
U
inself.source
the following is true:self.apply(U, mu).components(dofs) == restricted_op.apply(NumpyVectorArray(U.components(source_dofs)), mu))
Such an operator is mainly useful for
empirical interpolation
where the evaluation of the original operator only needs to be known for few selected degrees of freedom. If the operator has a small stencil, only fewsource_dofs
will be needed to evaluate the restricted operator which can make its evaluation very fast compared to evaluating the original operator.Parameters
- dofs
- One-dimensional
NumPy array
of degrees of freedom in the operatorrange
to which to restrict.
Returns
- restricted_op
- The restricted operator as defined above. The operator will have
NumpyVectorSpace
(len(source_dofs))
assource
andNumpyVectorSpace
(len(dofs))
asrange
. - source_dofs
- One-dimensional
NumPy array
of source degrees of freedom as defined above.
-
mpi module¶
-
class
pymor.operators.mpi.
MPIOperator
(obj_id, with_apply2=False, pickle_local_spaces=True, space_type=<class 'pymor.vectorarrays.mpi.MPIVectorSpace'>)[source]¶ Bases:
pymor.operators.basic.OperatorBase
MPI distributed
Operator
.Given a single-rank implementation of an
Operator
, this wrapper class uses the event loop frompymor.tools.mpi
to allow an MPI distributed usage of theOperator
.Instances of
MPIOperator
can be used on rank 0 like any other (non-distributed)Operator
.Note, however, that the underlying
Operator
implementation needs to be MPI aware. For instance, the operator’sapply
method has to perform the necessary MPI communication to obtain all DOFs hosted on other MPI ranks which are required for the local operator evaluation.Instead of instantiating
MPIOperator
directly, it is usually preferable to usempi_wrap_operator
instead.Parameters
- obj_id
ObjectId
of the localOperators
on each rank.- with_apply2
- Set to
True
if the operator implementation has its own MPI aware implementation ofapply2
andpairwise_apply2
. Otherwise, the default implementations usingapply
anddot
will be used. - pickle_local_spaces
- If
pickle_local_spaces
isFalse
, a unique identifier is computed for each local source/rangeVectorSpaceInterface
, which is then transferred to rank 0 instead of the trueVectorSpaceInterface
. This allows the useage ofMPIVectorArray
even when the localVectorSpaces
are not picklable. - space_type
- This class will be used to wrap the local
VectorArrays
returned by the local operators into an MPI distributedVectorArray
managed from rank 0. By default,MPIVectorSpace
will be used, other options areMPIVectorSpaceAutoComm
andMPIVectorSpaceNoComm
.
Methods
-
pymor.operators.mpi.
mpi_wrap_operator
(obj_id, with_apply2=False, pickle_local_spaces=True, space_type=<class 'pymor.vectorarrays.mpi.MPIVectorSpace'>)[source]¶ Wrap MPI distributed local
Operators
to a globalOperator
on rank 0.Given MPI distributed local
Operators
referred to by theObjectId
obj_id
, return a newOperator
which manages these distributed operators from rank 0. This is done by instantiatingMPIOperator
. Additionally, the structure of the wrapped operators is preserved. E.g.LincombOperators
will be wrapped as aLincombOperator
ofMPIOperators
.Parameters
See : class:
MPIOperator
.Returns
The wrapped
Operator
.
numpy module¶
This module provides the following NumPy
based Operators
:
NumpyMatrixOperator
wraps a 2DNumPy array
as anOperator
.NumpyMatrixBasedOperator
should be used as base class for allOperators
which assemble into aNumpyMatrixOperator
.NumpyGenericOperator
wraps an arbitrary Python function betweenNumPy arrays
as anOperator
.
-
class
pymor.operators.numpy.
NumpyGenericOperator
(mapping, transpose_mapping=None, dim_source=1, dim_range=1, linear=False, parameter_type=None, source_id=None, range_id=None, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.basic.OperatorBase
Wraps an arbitrary Python function between
NumPy arrays
as a anOperator
.Parameters
- mapping
The function to wrap. If
parameter_type
isNone
, the function is of the formmapping(U)
and is expected to be vectorized. In particular:mapping(U).shape == U.shape[:-1] + (dim_range,).
If
parameter_type
is notNone
, the function has to have the signaturemapping(U, mu)
.- transpose_mapping
The transpsoe function to wrap. If
parameter_type
isNone
, the function is of the formtranspose_mapping(U)
and is expected to be vectorized. In particular:transpose_mapping(U).shape == U.shape[:-1] + (dim_source,).
If
parameter_type
is notNone
, the function has to have the signaturetranspose_mapping(U, mu)
.- dim_source
- Dimension of the operator’s source.
- dim_range
- Dimension of the operator’s range.
- linear
- Set to
True
if the providedmapping
andtranspose_mapping
are linear. - parameter_type
- The
ParameterType
of theParameters
the mapping accepts. - name
- Name of the operator.
Methods
-
class
pymor.operators.numpy.
NumpyMatrixBasedOperator
[source]¶ Bases:
pymor.operators.basic.OperatorBase
Base class for operators which assemble into a
NumpyMatrixOperator
.Methods
Attributes
-
sparse
¶ True
if the operator assembles into a sparse matrix,False
if the operator assembles into a dense matrix,None
if unknown.
-
assemble
(mu=None)[source]¶ Assembles the operator for a given
Parameter
.Parameters
- mu
- The
Parameter
for which to assemble the operator.
Returns
The assembled parameter independent
Operator
.
-
export_matrix
(filename, matrix_name=None, output_format='matlab', mu=None)[source]¶ Save the matrix of the operator to a file.
Parameters
- filename
- Name of output file.
- matrix_name
- The name, the output matrix is given. (Comment field is used in
case of Matrix Market output_format.) If
None
, theOperator
‘sname
is used. - output_format
- Output file format. Either
matlab
ormatrixmarket
. - mu
- The
Parameter
to assemble the to be exported matrix for.
-
-
class
pymor.operators.numpy.
NumpyMatrixOperator
(matrix, source_id=None, range_id=None, solver_options=None, name=None)[source]¶ Bases:
pymor.operators.numpy.NumpyMatrixBasedOperator
Wraps a 2D
NumPy array
as anOperator
.Parameters
- matrix
- The
NumPy array
which is to be wrapped. - name
- Name of the operator.
Methods
Attributes
-
apply_inverse
(V, mu=None, least_squares=False, check_finite=True, default_sparse_solver_backend='scipy')[source]¶ Apply the inverse operator.
Parameters
- V
VectorArray
of vectors to which the inverse operator is applied.- mu
- The
Parameter
for which to evaluate the inverse operator. - least_squares
If
True
, solve the least squares problem:u = argmin ||op(u) - v||_2.
Since for an invertible operator the least squares solution agrees with the result of the application of the inverse operator, setting this option should, in general, have no effect on the result for those operators. However, note that when no appropriate
solver_options
are set for the operator, most implementations will choose a least squares solver by default which may be undesirable.- check_finite
- Test if solution only containes finite values.
- default_solver
- Default sparse solver backend to use (scipy, pyamg, generic).
Returns
VectorArray
of the inverse operator evaluations.Raises
- InversionError
- The operator could not be inverted.
Defaults
check_finite, default_sparse_solver_backend (see
pymor.core.defaults
)