Linear Operators¶
Module: pycsou.core.linop
Abstract classes for multidimensional linear maps.
|
Base class for linear operators. |
-
class
LinearOperator
(shape: Tuple[int, int], dtype: Optional[type] = None, is_explicit: bool = False, is_dense: bool = False, is_sparse: bool = False, is_dask: bool = False, is_symmetric: bool = False, lipschitz_cst: float = inf)[source]¶ Bases:
pycsou.core.map.DifferentiableMap
Base class for linear operators.
Any instance/subclass of this class must at least implement the abstract methods
__call__
andadjoint
.Notes
This class supports the following arithmetic operators
+
,-
,*
,@
,**
and/
, implemented with the class methods__add__
/__radd__
,__sub__
/__neg__
,__mul__
/__rmul__
,__matmul__
,__pow__
,__truediv__
. Such arithmetic operators can be used to add, substract, scale, compose, exponentiate or evaluateLinearOperator
instances.-
__init__
(shape: Tuple[int, int], dtype: Optional[type] = None, is_explicit: bool = False, is_dense: bool = False, is_sparse: bool = False, is_dask: bool = False, is_symmetric: bool = False, lipschitz_cst: float = inf)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
matvec
(x: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Alias for
self.__call__
to comply with Scipy’s interface.
-
abstract
adjoint
(y: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-
transpose
(y: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the tranpose of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the transpose should be evaluated.
- Returns
Evaluation of the transpose at
y
.- Return type
Union[Number, np.ndarray]
Notes
For real-valued operators, the adjoint and the transpose coincide. For complex-valued operators, we can define the transpose as \(\mathbf{A}^T \mathbf{y}=\overline{\mathbf{A}^\ast \overline{\mathbf{y}}}\).
-
jacobianT
(arg: Union[numbers.Number, numpy.ndarray, None] = None) → pycsou.core.linop.LinearOperator[source]¶ Transpose of the Jacobian matrix of the differentiable map evaluated at
arg
.- Parameters
arg (Union[Number, np.ndarray]) – Point at which the transposed Jacobian matrix is evaluated.
- Returns
Linear operator associated to the transposed Jacobian matrix.
- Return type
-
get_adjointOp
() → Union[pycsou.core.linop.LinearOperator, pycsou.core.linop.AdjointLinearOperator][source]¶ Return the adjoint operator as a
LinearOperator
instance.- Returns
The adjoint operator.
- Return type
Union[‘LinearOperator’, ‘AdjointLinearOperator’]
-
property
H
¶ Alias for
self.get_adjointOp
.
-
get_transposeOp
() → pycsou.core.linop.TransposeLinearOperator[source]¶ Return the transpose operator as a
LinearOperator
instance.- Returns
The transpose operator.
- Return type
Notes
For real-valued operators, the adjoint and the transpose coincide. For complex-valued operators, we can define the transpose as \(\mathbf{A}^T \mathbf{y}=\overline{\mathbf{A}^\ast \overline{\mathbf{y}}}\).
-
property
T
¶ Alias for
self.get_adjointOp
.
-
property
RangeGram
¶ Range-Gram operator, obtained by the composition
self * self.H
(adjoint followed by operator).- Returns
The Range-Gram operator.
- Return type
-
property
DomainGram
¶ Domain-Gram operator, obtained by the composition
self.H * self
(operator followed by adjoint).- Returns
The Domain-Gram operator.
- Return type
-
eigenvals
(k: int, which='LM', **kwargs) → numpy.ndarray[source]¶ Find
k
eigenvalues of a square operator.- Parameters
k (int) – The number of eigenvalues and eigenvectors desired.
k
must be strictly smaller than the dimension of the operator. It is not possible to compute all eigenvectors of a matrix.which (str, [‘LM’ | ‘SM’ | ‘LR’ | ‘SR’ | ‘LI’ | ‘SI’]) –
Which
k
eigenvalues to find:‘LM’ : largest magnitude
‘SM’ : smallest magnitude
‘LR’ : largest real part
‘SR’ : smallest real part
‘LI’ : largest imaginary part
‘SI’ : smallest imaginary part
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s functions
scipy.sparse.linalg.eigs()
andscipy.sparse.linalg.eigsh()
.
- Returns
Array containing the
k
requested eigenvalues.- Return type
np.ndarray
- Raises
NotImplementedError – If the linear operator is not square.
Notes
This function calls one of the two Scipy’s functions:
scipy.sparse.linalg.eigs()
andscipy.sparse.linalg.eigsh()
. See the documentation of these two functions for more information on their behaviour and the underlying ARPACK functions they rely on.See also
svds()
-
singularvals
(k: int, which='LM', **kwargs) → numpy.ndarray[source]¶ Compute the largest or smallest
k
singular values of an operator. The order of the singular values is not guaranteed.- Parameters
k (int) – Number of singular values to compute. Must be
1 <= k < min(self.shape)
.which (str, [‘LM’ | ‘SM’]) –
Which
k
singular values to find:‘LM’ : largest magnitude
‘SM’ : smallest magnitude
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s function
scipy.sparse.linalg.svds()
.
- Returns
Array containing the
k
requested singular values.- Return type
np.ndarray
Examples
>>> sig = np.repeat([0., 1., 0.], 10) >>> filter = signal.hann(5); filter[filter.size//2:] = 0 >>> ConvOp = Convolve1D(size=sig.size, filter=filter) >>> np.round((ConvOp.singularvals(k=3, which='LM', tol=1e-3)), 2) array([0.5, 0.5, 0.5])
Notes
This function calls the Scipy’s function:
scipy.sparse.linalg.svds()
. See the documentation of this function for more information on its behaviour and the underlying ARPACK/LOBPCG functions it relies on.See also
-
compute_lipschitz_cst
(**kwargs)[source]¶ Compute the Lipschitz constant of the operator.
- Parameters
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s functions
scipy.sparse.linalg.eigs()
,scipy.sparse.linalg.eigsh()
,scipy.sparse.linalg.svds()
.- Returns
Nothing: The Lipschitz constant is stored in the attribute
self.lipschitz_cst
.- Return type
Examples
>>> sig = np.repeat([0., 1., 0.], 10) >>> filter = signal.hann(5); filter[filter.size//2:] = 0 >>> ConvOp = Convolve1D(size=sig.size, filter=filter) >>> ConvOp.compute_lipschitz_cst(tol=1e-2); np.round(ConvOp.lipschitz_cst,1) 0.5
Notes
The Lipschtiz constant of a linear operator is its largest singular value. This function therefore calls the methods
self.singularvals
withk=1
andwhich='LM'
to perform this computation.Warning
For high-dimensional linear operators this method can be very time-consuming. Reducing the computation accuracy with the optional argument
tol: float
may help reduce the computational burden. See Scipy’s functionsscipy.sparse.linalg.eigs()
,scipy.sparse.linalg.eigsh()
,scipy.sparse.linalg.svds()
for more on this parameter.
-
todense
() → DenseLinearOperator[source]¶ Convert the operator to a
DenseLinearOperator
.- Returns
The dense linear operator representation.
- Return type
-
tosparse
() → SparseLinearOperator[source]¶ Convert the operator to a
SparseLinearOperator
.- Returns
The sparse linear operator representation.
- Return type
-
tosciop
() → scipy.sparse.linalg.interface.LinearOperator[source]¶ Convert the operator to a Scipy
scipy.sparse.linalg.LinearOperator
.- Returns
The Scipy linear operator representation.
- Return type
-
property
SciOp
¶ Alias for method
self.tosciop
.
-
topylop
() → pylops.LinearOperator.LinearOperator[source]¶ Convert the operator to a Pylops
pylops.LinearOperator
.- Returns
The Pylops linear operator representation.
- Return type
-
property
PyLop
¶ Alias for method
self.topylop
.
-
cond
(**kwargs) → float[source]¶ Compute the condition number of the operator.
- Parameters
kwargs (dict) – A dict of additional keyword arguments values accepted by Pylops’ method
pylops.LinearOperator.cond()
.- Returns
Condition number.
- Return type
-
pinv
(y: Union[numbers.Number, numpy.ndarray], eps: numbers.Number = 0, **kwargs) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluate the pseudo-inverse of the operator at
y
.- Parameters
y (Union[Number, np.ndarray]) – Point at which the pseudo-inverse is evaluated.
eps (Number) – Tikhonov damping.
kwargs – Arbitrary keyword arguments accepted by the function:
pylops.optimization.leastsquares.NormalEquationsInversion()
.
- Returns
Evaluation of the pseudo-inverse of the operator at
y
.- Return type
Notes
This is a wrapper around the function
pylops.optimization.leastsquares.NormalEquationsInversion()
. Additional information can be found in the help of this function.
-
property
PinvOp
¶ Return the pseudo-inverse of the operator as a
LinearOperator
.
-
property
dagger
¶ Alias for
self.PinvOp
.
-
property
RowProjector
¶ Orthogonal projection operator onto the rows of
self
. It is given byself.dagger * self
.
-
property
ColProjector
¶ Orthogonal projection operator onto the columns of
self
. It is given byself * self.dagger
.
-
-
class
AdjointLinearOperator
(LinOp: pycsou.core.linop.LinearOperator)[source]¶ Bases:
pycsou.core.linop.LinearOperator
-
__init__
(LinOp: pycsou.core.linop.LinearOperator)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
adjoint
(x: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-
compute_lipschitz_cst
(**kwargs)[source]¶ Compute the Lipschitz constant of the operator.
- Parameters
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s functions
scipy.sparse.linalg.eigs()
,scipy.sparse.linalg.eigsh()
,scipy.sparse.linalg.svds()
.- Returns
Nothing: The Lipschitz constant is stored in the attribute
self.lipschitz_cst
.- Return type
Examples
>>> sig = np.repeat([0., 1., 0.], 10) >>> filter = signal.hann(5); filter[filter.size//2:] = 0 >>> ConvOp = Convolve1D(size=sig.size, filter=filter) >>> ConvOp.compute_lipschitz_cst(tol=1e-2); np.round(ConvOp.lipschitz_cst,1) 0.5
Notes
The Lipschtiz constant of a linear operator is its largest singular value. This function therefore calls the methods
self.singularvals
withk=1
andwhich='LM'
to perform this computation.Warning
For high-dimensional linear operators this method can be very time-consuming. Reducing the computation accuracy with the optional argument
tol: float
may help reduce the computational burden. See Scipy’s functionsscipy.sparse.linalg.eigs()
,scipy.sparse.linalg.eigsh()
,scipy.sparse.linalg.svds()
for more on this parameter.
-
-
class
TransposeLinearOperator
(LinOp: pycsou.core.linop.LinearOperator)[source]¶ Bases:
pycsou.core.linop.LinearOperator
-
__init__
(LinOp: pycsou.core.linop.LinearOperator)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
adjoint
(x: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-
-
class
LinOpSum
(LinOp1: pycsou.core.linop.LinearOperator, LinOp2: pycsou.core.linop.LinearOperator, dtype: Optional[type] = None)[source]¶ Bases:
pycsou.core.linop.LinearOperator
,pycsou.core.map.DiffMapSum
-
__init__
(LinOp1: pycsou.core.linop.LinearOperator, LinOp2: pycsou.core.linop.LinearOperator, dtype: Optional[type] = None)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
adjoint
(y: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-
-
class
LinOpComp
(LinOp1: pycsou.core.linop.LinearOperator, LinOp2: pycsou.core.linop.LinearOperator, dtype: Optional[type] = None)[source]¶ Bases:
pycsou.core.linop.LinearOperator
,pycsou.core.map.DiffMapComp
-
__init__
(LinOp1: pycsou.core.linop.LinearOperator, LinOp2: pycsou.core.linop.LinearOperator, dtype: Optional[type] = None)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
adjoint
(y: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-
-
class
SymmetricLinearOperator
(LinOp: pycsou.core.linop.LinearOperator)[source]¶ Bases:
pycsou.core.linop.LinearOperator
-
__init__
(LinOp: pycsou.core.linop.LinearOperator)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
adjoint
(y: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-
-
class
UnitaryOperator
(size: int, dtype: Optional[type] = None, is_explicit: bool = False, is_dense: bool = False, is_sparse: bool = False, is_dask: bool = False, is_symmetric: bool = False)[source]¶ Bases:
pycsou.core.linop.LinearOperator
-
__init__
(size: int, dtype: Optional[type] = None, is_explicit: bool = False, is_dense: bool = False, is_sparse: bool = False, is_dask: bool = False, is_symmetric: bool = False)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
property
RangeGram
¶ Range-Gram operator, obtained by the composition
self * self.H
(adjoint followed by operator).- Returns
The Range-Gram operator.
- Return type
-
property
DomainGram
¶ Domain-Gram operator, obtained by the composition
self.H * self
(operator followed by adjoint).- Returns
The Domain-Gram operator.
- Return type
-
eigenvals
(k: int, which='LM', **kwargs) → numpy.ndarray[source]¶ Find
k
eigenvalues of a square operator.- Parameters
k (int) – The number of eigenvalues and eigenvectors desired.
k
must be strictly smaller than the dimension of the operator. It is not possible to compute all eigenvectors of a matrix.which (str, [‘LM’ | ‘SM’ | ‘LR’ | ‘SR’ | ‘LI’ | ‘SI’]) –
Which
k
eigenvalues to find:‘LM’ : largest magnitude
‘SM’ : smallest magnitude
‘LR’ : largest real part
‘SR’ : smallest real part
‘LI’ : largest imaginary part
‘SI’ : smallest imaginary part
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s functions
scipy.sparse.linalg.eigs()
andscipy.sparse.linalg.eigsh()
.
- Returns
Array containing the
k
requested eigenvalues.- Return type
np.ndarray
- Raises
NotImplementedError – If the linear operator is not square.
Notes
This function calls one of the two Scipy’s functions:
scipy.sparse.linalg.eigs()
andscipy.sparse.linalg.eigsh()
. See the documentation of these two functions for more information on their behaviour and the underlying ARPACK functions they rely on.See also
svds()
-
singularvals
(k: int, which='LM', **kwargs) → numpy.ndarray[source]¶ Compute the largest or smallest
k
singular values of an operator. The order of the singular values is not guaranteed.- Parameters
k (int) – Number of singular values to compute. Must be
1 <= k < min(self.shape)
.which (str, [‘LM’ | ‘SM’]) –
Which
k
singular values to find:‘LM’ : largest magnitude
‘SM’ : smallest magnitude
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s function
scipy.sparse.linalg.svds()
.
- Returns
Array containing the
k
requested singular values.- Return type
np.ndarray
Examples
>>> sig = np.repeat([0., 1., 0.], 10) >>> filter = signal.hann(5); filter[filter.size//2:] = 0 >>> ConvOp = Convolve1D(size=sig.size, filter=filter) >>> np.round((ConvOp.singularvals(k=3, which='LM', tol=1e-3)), 2) array([0.5, 0.5, 0.5])
Notes
This function calls the Scipy’s function:
scipy.sparse.linalg.svds()
. See the documentation of this function for more information on its behaviour and the underlying ARPACK/LOBPCG functions it relies on.See also
-
compute_lipschitz_cst
(**kwargs)[source]¶ Compute the Lipschitz constant of the operator.
- Parameters
kwargs (dict) – A dict of additional keyword arguments values accepted by Scipy’s functions
scipy.sparse.linalg.eigs()
,scipy.sparse.linalg.eigsh()
,scipy.sparse.linalg.svds()
.- Returns
Nothing: The Lipschitz constant is stored in the attribute
self.lipschitz_cst
.- Return type
Examples
>>> sig = np.repeat([0., 1., 0.], 10) >>> filter = signal.hann(5); filter[filter.size//2:] = 0 >>> ConvOp = Convolve1D(size=sig.size, filter=filter) >>> ConvOp.compute_lipschitz_cst(tol=1e-2); np.round(ConvOp.lipschitz_cst,1) 0.5
Notes
The Lipschtiz constant of a linear operator is its largest singular value. This function therefore calls the methods
self.singularvals
withk=1
andwhich='LM'
to perform this computation.Warning
For high-dimensional linear operators this method can be very time-consuming. Reducing the computation accuracy with the optional argument
tol: float
may help reduce the computational burden. See Scipy’s functionsscipy.sparse.linalg.eigs()
,scipy.sparse.linalg.eigsh()
,scipy.sparse.linalg.svds()
for more on this parameter.
-
pinv
(y: Union[numbers.Number, numpy.ndarray], eps: numbers.Number = 0, **kwargs) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluate the pseudo-inverse of the operator at
y
.- Parameters
y (Union[Number, np.ndarray]) – Point at which the pseudo-inverse is evaluated.
eps (Number) – Tikhonov damping.
kwargs – Arbitrary keyword arguments accepted by the function:
pylops.optimization.leastsquares.NormalEquationsInversion()
.
- Returns
Evaluation of the pseudo-inverse of the operator at
y
.- Return type
Notes
This is a wrapper around the function
pylops.optimization.leastsquares.NormalEquationsInversion()
. Additional information can be found in the help of this function.
-
property
PinvOp
¶ Return the pseudo-inverse of the operator as a
LinearOperator
.
-
cond
(**kwargs)[source]¶ Compute the condition number of the operator.
- Parameters
kwargs (dict) – A dict of additional keyword arguments values accepted by Pylops’ method
pylops.LinearOperator.cond()
.- Returns
Condition number.
- Return type
-
-
class
LinOpPinv
(LinOp: pycsou.core.linop.LinearOperator, eps: numbers.Number = 0)[source]¶ Bases:
pycsou.core.linop.LinearOperator
-
__init__
(LinOp: pycsou.core.linop.LinearOperator, eps: numbers.Number = 0)[source]¶ - Parameters
dtype (Optional[type]) – Data type of the linear operator.
is_explicit (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy/Scipy/Dask array.is_dense (bool) – If
True
, the linear operator is specified explicitly in terms of a Numpy array.is_sparse (bool) – If
True
, the linear operator is specified explicitly in terms of a Scipy sparse matrix.is_dask (bool) – If
True
, the linear operator is specified explicitly in terms of a Dask array.is_symmetric (bool) – Whether the linear operator is symmetric or not.
lipschitz_cst (float) – Lipschitz constant (maximal singular value) of the linear operator if known. Default to \(+\infty\).
-
adjoint
(y: Union[numbers.Number, numpy.ndarray]) → Union[numbers.Number, numpy.ndarray][source]¶ Evaluates the adjoint of the operator at a point.
- Parameters
y (Union[Number, np.ndarray]) – Point at which the adjoint should be evaluated.
- Returns
Evaluation of the adjoint at
y
.- Return type
Union[Number, np.ndarray]
-