pymor.operators.numpy¶
Operators based on NumPy arrays.
This module provides the following NumPy-based Operators:
NumpyMatrixOperatorwraps a 2DNumPy arrayas anOperator.NumpyMatrixBasedOperatorshould be used as base class for allOperatorswhich assemble into aNumpyMatrixOperator.NumpyGenericOperatorwraps an arbitrary Python function betweenNumPy arraysas anOperator.
Module Contents¶
Classes¶
Wraps an arbitrary Python function between |
|
Base class for operators which assemble into a |
|
Wraps a 2D |
- class pymor.operators.numpy.NumpyGenericOperator(mapping, adjoint_mapping=None, dim_source=1, dim_range=1, linear=False, parameters={}, source_id=None, range_id=None, solver_options=None, name=None)[source]¶
Bases:
pymor.operators.interface.OperatorWraps an arbitrary Python function between
NumPy arraysas anOperator.Parameters
- mapping
The function to wrap. If
parametersisNone, the function is of the formmapping(U)and is expected to be vectorized. In particular:mapping(U).shape == U.shape[:-1] + (dim_range,).
If
parametersis notNone, the function has to have the signaturemapping(U, mu).- adjoint_mapping
The adjoint function to wrap. If
parametersisNone, the function is of the formadjoint_mapping(U)and is expected to be vectorized. In particular:adjoint_mapping(U).shape == U.shape[:-1] + (dim_source,).
If
parametersis notNone, the function has to have the signatureadjoint_mapping(U, mu).- dim_source
Dimension of the operator’s source.
- dim_range
Dimension of the operator’s range.
- linear
Set to
Trueif the providedmappingandadjoint_mappingare linear.- parameters
The
Parametersthe depends on.- solver_options
The
solver_optionsfor the operator.- name
Name of the operator.
- apply(self, U, mu=None)[source]¶
Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
- apply_adjoint(self, V, mu=None)[source]¶
Apply the adjoint operator.
For any given linear
Operatorop,parameter valuesmuandVectorArraysU,Vin thesourceresp.rangewe have:op.apply_adjoint(V, mu).dot(U) == V.inner(op.apply(U, mu))
Thus, when
opis represented by a matrixM,apply_adjointis given by left-multplication of (the complex conjugate of)MwithV.Parameters
- V
VectorArrayof vectors to which the adjoint operator is applied.- mu
The
parameter valuesfor which to apply the adjoint operator.
Returns
VectorArrayof the adjoint operator evaluations.
- class pymor.operators.numpy.NumpyMatrixBasedOperator[source]¶
Bases:
pymor.operators.interface.OperatorBase class for operators which assemble into a
NumpyMatrixOperator.- sparse[source]¶
Trueif the operator assembles into a sparse matrix,Falseif the operator assembles into a dense matrix,Noneif unknown.
- assemble(self, mu=None)[source]¶
Assemble the operator for given
parameter values.The result of the method strongly depends on the given operator. For instance, a matrix-based operator will assemble its matrix, a
LincombOperatorwill 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 valuesfor which to assemble the operator.
Returns
Parameter-independent, assembled
Operator.
- apply(self, U, mu=None)[source]¶
Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
- apply_adjoint(self, V, mu=None)[source]¶
Apply the adjoint operator.
For any given linear
Operatorop,parameter valuesmuandVectorArraysU,Vin thesourceresp.rangewe have:op.apply_adjoint(V, mu).dot(U) == V.inner(op.apply(U, mu))
Thus, when
opis represented by a matrixM,apply_adjointis given by left-multplication of (the complex conjugate of)MwithV.Parameters
- V
VectorArrayof vectors to which the adjoint operator is applied.- mu
The
parameter valuesfor which to apply the adjoint operator.
Returns
VectorArrayof the adjoint operator evaluations.
- as_range_array(self, mu=None)[source]¶
Return a
VectorArrayrepresentation of the operator in its range space.In the case of a linear operator with
NumpyVectorSpaceassource, this method returns for givenparameter valuesmuaVectorArrayVin the operator’srange, such thatV.lincomb(U.to_numpy()) == self.apply(U, mu)
for all
VectorArraysU.Parameters
- mu
The
parameter valuesfor which to return theVectorArrayrepresentation.
Returns
- V
The
VectorArraydefined above.
- as_source_array(self, mu=None)[source]¶
Return a
VectorArrayrepresentation of the operator in its source space.In the case of a linear operator with
NumpyVectorSpaceasrange, this method returns for givenparameter valuesmuaVectorArrayVin the operator’ssource, such thatself.range.make_array(V.inner(U).T) == self.apply(U, mu)
for all
VectorArraysU.Parameters
- mu
The
parameter valuesfor which to return theVectorArrayrepresentation.
Returns
- V
The
VectorArraydefined above.
- apply_inverse(self, V, mu=None, initial_guess=None, least_squares=False)[source]¶
Apply the inverse operator.
Parameters
- V
VectorArrayof vectors to which the inverse operator is applied.- mu
The
parameter valuesfor which to evaluate the inverse operator.- initial_guess
VectorArraywith the same length asVcontaining initial guesses for the solution. Some implementations ofapply_inversemay ignore this parameter. IfNonea solver-dependent default is used.- 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_optionsare set for the operator, most implementations will choose a least squares solver by default which may be undesirable.
Returns
VectorArrayof the inverse operator evaluations.Raises
- InversionError
The operator could not be inverted.
- export_matrix(self, 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’snameis used.- output_format
Output file format. Either
matlabormatrixmarket.- mu
The
parameter valuesto 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:
NumpyMatrixBasedOperatorWraps a 2D
NumPy arrayas anOperator.Parameters
- matrix
The
NumPy arraywhich is to be wrapped.- source_id
The id of the operator’s
sourceVectorSpace.- range_id
The id of the operator’s
rangeVectorSpace.- solver_options
The
solver_optionsfor the operator.- name
Name of the operator.
- classmethod from_file(cls, path, key=None, source_id=None, range_id=None, solver_options=None, name=None)[source]¶
- assemble(self, mu=None)[source]¶
Assemble the operator for given
parameter values.The result of the method strongly depends on the given operator. For instance, a matrix-based operator will assemble its matrix, a
LincombOperatorwill 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 valuesfor which to assemble the operator.
Returns
Parameter-independent, assembled
Operator.
- as_range_array(self, mu=None)[source]¶
Return a
VectorArrayrepresentation of the operator in its range space.In the case of a linear operator with
NumpyVectorSpaceassource, this method returns for givenparameter valuesmuaVectorArrayVin the operator’srange, such thatV.lincomb(U.to_numpy()) == self.apply(U, mu)
for all
VectorArraysU.Parameters
- mu
The
parameter valuesfor which to return theVectorArrayrepresentation.
Returns
- V
The
VectorArraydefined above.
- as_source_array(self, mu=None)[source]¶
Return a
VectorArrayrepresentation of the operator in its source space.In the case of a linear operator with
NumpyVectorSpaceasrange, this method returns for givenparameter valuesmuaVectorArrayVin the operator’ssource, such thatself.range.make_array(V.inner(U).T) == self.apply(U, mu)
for all
VectorArraysU.Parameters
- mu
The
parameter valuesfor which to return theVectorArrayrepresentation.
Returns
- V
The
VectorArraydefined above.
- apply(self, U, mu=None)[source]¶
Apply the operator to a
VectorArray.Parameters
- U
VectorArrayof vectors to which the operator is applied.- mu
The
parameter valuesfor which to evaluate the operator.
Returns
VectorArrayof the operator evaluations.
- apply_adjoint(self, V, mu=None)[source]¶
Apply the adjoint operator.
For any given linear
Operatorop,parameter valuesmuandVectorArraysU,Vin thesourceresp.rangewe have:op.apply_adjoint(V, mu).dot(U) == V.inner(op.apply(U, mu))
Thus, when
opis represented by a matrixM,apply_adjointis given by left-multplication of (the complex conjugate of)MwithV.Parameters
- V
VectorArrayof vectors to which the adjoint operator is applied.- mu
The
parameter valuesfor which to apply the adjoint operator.
Returns
VectorArrayof the adjoint operator evaluations.
- apply_inverse(self, V, mu=None, initial_guess=None, least_squares=False, check_finite=True, default_sparse_solver_backend='scipy')[source]¶
Apply the inverse operator.
Parameters
- V
VectorArrayof vectors to which the inverse operator is applied.- mu
The
parameter valuesfor which to evaluate the inverse operator.- initial_guess
VectorArraywith the same length asVcontaining initial guesses for the solution. Some implementations ofapply_inversemay ignore this parameter. IfNonea solver-dependent default is used.- 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_optionsare set for the operator, most implementations will choose a least squares solver by default which may be undesirable.- check_finite
Test if solution only contains finite values.
- default_sparse_solver_backend
Default sparse solver backend to use (scipy, pyamg, generic).
Returns
VectorArrayof the inverse operator evaluations.Raises
- InversionError
The operator could not be inverted.
- apply_inverse_adjoint(self, U, mu=None, initial_guess=None, least_squares=False)[source]¶
Apply the inverse adjoint operator.
Parameters
- U
VectorArrayof vectors to which the inverse adjoint operator is applied.- mu
The
parameter valuesfor which to evaluate the inverse adjoint operator.- initial_guess
VectorArraywith the same length asUcontaining initial guesses for the solution. Some implementations ofapply_inverse_adjointmay ignore this parameter. IfNonea solver-dependent default is used.- 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_optionsare set for the operator, most operator implementations will choose a least squares solver by default which may be undesirable.
Returns
VectorArrayof the inverse adjoint operator evaluations.Raises
- InversionError
The operator could not be inverted.
- _assemble_lincomb(self, operators, coefficients, identity_shift=0.0, solver_options=None, name=None)[source]¶
Try to assemble a linear combination of the given operators.
Returns a new
Operatorwhich represents the sumc_1*O_1 + ... + c_N*O_N + s*I
where
O_iareOperators,c_i,sscalar coefficients andIthe identity.This method is called in the
assemblemethod ofLincombOperatoron the first of its operators. 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 returnsNoneto indicate that assembly is not possible.Parameters
- operators
List of
OperatorsO_iwhose linear combination is formed.- coefficients
List of the corresponding linear coefficients
c_i.- identity_shift
The coefficient
s.- solver_options
solver_optionsfor the assembled operator.- name
Name of the assembled operator.
Returns
The assembled
Operatorif assembly is possible, otherwiseNone.