linalg.basic#

provides the fundamental classes for linear algebra operations.

Classes#

Inheritance diagram of psydac.linalg.basic

ComposedLinearOperator(domain, codomain, *args)

Composition \(A_n\circ\dots\circ A_1\) of two or more linear operators \(A_1,\dots,A_n\).

IdentityOperator(domain[, codomain])

Identity operator acting between a vector space V and itself.

InverseLinearOperator(A, **kwargs)

Abstract base class for the (approximate) inverse \(A^{-1}\) of a square matrix \(A\).

LinearOperator()

Abstract base class for all linear operators acting between two vector spaces V (domain) and W (codomain).

LinearSolver()

Solver for the square linear system Ax=b, where x and b belong to the same vector space V.

MatrixFreeLinearOperator(domain, codomain, dot)

General linear operator represented by a callable dot method.

PowerLinearOperator(domain, codomain, A, n)

Power \(A^n\) of a linear operator \(A\) for some integer \(n\geq 0\).

ScaledLinearOperator(domain, codomain, c, A)

A linear operator \(A\) scalar multiplied by a real constant \(c\).

SumLinearOperator(domain, codomain, *args)

Sum \(\sum_{i=1}^n A_i\) of linear operators \(A_1,\dots,A_n\) acting between the same vector spaces V (domain) and W (codomain).

Vector()

Element of a vector space V.

VectorSpace()

Finite-dimensional vector space V with a scalar (dot) product.

ZeroOperator(domain[, codomain])

Zero operator mapping any vector from its domain V to the zero vector of its codomain W.

Details#

provides the fundamental classes for linear algebra operations.

class VectorSpace[source]#

Bases: ABC

Finite-dimensional vector space V with a scalar (dot) product.

abstract property dimension#

The dimension of a vector space V is the cardinality (i.e. the number of vectors) of a basis of V over its base field.

abstract property dtype#

The data type of the field over which the space is built.

abstract zeros()[source]#

Get a copy of the null element of the vector space V.

Returns:
nullVector

A new vector object with all components equal to zero.

dot(a, b)[source]#

Evaluate the scalar product between two vectors of the same space.

abstract axpy(a, x, y)[source]#

Increment the vector y with the a-scaled vector x, i.e. y = a * x + y, provided that x and y belong to the same vector space V (self). The scalar value a may be real or complex, depending on the field of V.

Parameters:
ascalar

The scaling coefficient needed for the operation.

xVector

The vector which is not modified by this function.

yVector

The vector modified by this function (incremented by a * x).

class Vector[source]#

Bases: ABC

Element of a vector space V.

property shape#

A tuple containing the dimension of the space.

property dtype#

The data type of the vector field V this vector belongs to.

dot(v)[source]#

Evaluate the scalar product with the vector v of the same space.

Parameters:
vVector

Vector belonging to the same space as self.

mul_iadd(a, v)[source]#

Compute self += a * v, where v is another vector of the same space.

Parameters:
ascalar

Rescaling coefficient, which can be cast to the correct dtype.

vVector

Vector belonging to the same space as self.

abstract property space#

Vector space to which this vector belongs.

abstract toarray(**kwargs)[source]#

Convert to Numpy 1D array.

abstract copy(out=None)[source]#

Ensure x.copy(out=x) returns x and not a new object.

abstract conjugate(out=None)[source]#

Compute the complex conjugate vector.

If the field is real (i.e. self.dtype in (np.float32, np.float64)) this method is equivalent to copy. If the field is complex (i.e. self.dtype in (np.complex64, np.complex128)) this method returns the complex conjugate of self, element-wise.

The behavior of this function is similar to numpy.conjugate(self, out=None).

conj(out=None)[source]#

Compute the complex conjugate vector.

If the field is real (i.e. self.dtype in (np.float32, np.float64)) this method is equivalent to copy. If the field is complex (i.e. self.dtype in (np.complex64, np.complex128)) this method returns the complex conjugate of self, element-wise.

The behavior of this function is similar to numpy.conj(self, out=None).

class LinearOperator[source]#

Bases: ABC

Abstract base class for all linear operators acting between two vector spaces V (domain) and W (codomain).

property shape#

A tuple containing the dimension of the codomain and domain.

abstract property domain#

The domain of the linear operator - an element of Vectorspace

abstract property codomain#

The codomain of the linear operator - an element of Vectorspace

abstract property dtype#
abstract tosparse()[source]#
abstract toarray()[source]#

Convert to Numpy 2D array.

abstract dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

abstract transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

property T#

Calls transpose method to return the transpose of self.

property H#

Calls transpose method with conjugate=True flag to return the Hermitian transpose of self.

idot(v, out)[source]#

Implements out += self @ v with a temporary. Subclasses should provide an implementation without a temporary.

class ZeroOperator(domain, codomain=None)[source]#

Bases: LinearOperator

Zero operator mapping any vector from its domain V to the zero vector of its codomain W.

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property dtype#
copy()[source]#
toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

class IdentityOperator(domain, codomain=None)[source]#

Bases: LinearOperator

Identity operator acting between a vector space V and itself. Useful for example in custom linear operator classes together with the apply_essential_bc method to create projection operators.

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property dtype#
copy()[source]#

Returns a new IdentityOperator object acting between the same vector spaces.

toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Could return self, but by convention returns new object.

dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

class ScaledLinearOperator(domain, codomain, c, A)[source]#

Bases: LinearOperator

A linear operator \(A\) scalar multiplied by a real constant \(c\).

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property scalar#

Returns the scalar value by which the operator is multiplied.

property operator#

Returns the operator that is multiplied by the scalar.

property dtype#
toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

class SumLinearOperator(domain, codomain, *args)[source]#

Bases: LinearOperator

Sum \(\sum_{i=1}^n A_i\) of linear operators \(A_1,\dots,A_n\) acting between the same vector spaces V (domain) and W (codomain).

property domain#

The domain of the linear operator, element of class VectorSpace.

property codomain#

The codomain of the linear operator, element of class VectorSpace.

property addends#

A tuple containing the addends of the linear operator, elements of class LinearOperator.

property dtype#
toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

static simplify(addends)[source]#

Simplifies a sum of linear operators by combining addends of the same class.

dot(v, out=None)[source]#

Evaluates SumLinearOperator object at a vector v element of domain.

class ComposedLinearOperator(domain, codomain, *args)[source]#

Bases: LinearOperator

Composition \(A_n\circ\dots\circ A_1\) of two or more linear operators \(A_1,\dots,A_n\).

property tmp_vectors#

A tuple containing the storage vectors that are repeatedly being used upon calling the dot method. This avoids the creation of new vectors at each call of the dot method.

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property multiplicants#

A tuple \((A_1,\dots,A_n)\) containing the multiplicants of the linear operator \(self = A_n\circ\dots\circ A_1\).

property dtype#
toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

exchange_assembly_data()[source]#
set_backend(backend)[source]#
class PowerLinearOperator(domain, codomain, A, n)[source]#

Bases: LinearOperator

Power \(A^n\) of a linear operator \(A\) for some integer \(n\geq 0\).

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property dtype#
property operator#

Returns the operator that is raised to the power.

property factorial#

Returns the power to which the operator is raised.

toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

class InverseLinearOperator(A, **kwargs)[source]#

Bases: LinearOperator

Abstract base class for the (approximate) inverse \(A^{-1}\) of a square matrix \(A\). The result of A_inv.dot(b) is the (approximate) solution x of the linear system A x = b, where x and b belong to the same vector space V.

We assume that the linear system is solved by an iterative method, which needs a first guess x0 and an exit condition based on tol and maxiter.

Concrete subclasses of this class must implement the dot method and take care of any internal storage which might be necessary.

Parameters:
Apsydac.linalg.basic.LinearOperator

Left-hand-side matrix A of linear system.

x0psydac.linalg.basic.Vector

First guess of solution for iterative solver (optional).

tolfloat

Absolute tolerance for L2-norm of residual r = A*x - b.

maxiter: int

Maximum number of iterations.

verbosebool

If True, L2-norm of residual r is printed at each iteration.

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property dtype#
property linop#

The linear operator \(A\) of which this object is the inverse \(A^{-1}\).

The linear operator \(A\) can be modified in place, or replaced entirely through the setter. A substitution should only be made in cases where no other options are viable, as it breaks the one-to-one map between the original linear operator \(A\) (passed to the constructor) and the current InverseLinearOperator object \(A^{-1}\). Use with extreme care!

toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
get_info()[source]#

Returns the previous convergence information.

get_options(key=None)[source]#

Get a copy of all the solver options, or a specific value of interest.

Parameters:
keystr | None

Name of the specific option of interest (default: None).

Returns:
dict | type(self._options[‘key’]) | None

If key is given, get the specific option of interest. If there is no such option, None is returned instead. If key is not given, get a copy of all the solver options in a dictionary.

set_options(**kwargs)[source]#

Set the solver options by passing keyword arguments.

transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

class LinearSolver[source]#

Bases: ABC

Solver for the square linear system Ax=b, where x and b belong to the same vector space V.

property shape#
abstract property space#
abstract transpose()[source]#

Return the transpose of the LinearSolver.

abstract solve(rhs, out=None)[source]#
property T#
class MatrixFreeLinearOperator(domain, codomain, dot, dot_transpose=None)[source]#

Bases: LinearOperator

General linear operator represented by a callable dot method.

Parameters:
domainVectorSpace

The domain of the linear operator.

codomainVectorSpace

The codomain of the linear operator.

dotCallable

The method of the linear operator, assumed to map from domain to codomain. This method can take out as an optional argument but this is not mandatory.

dot_transpose: Callable

The method of the transpose of the linear operator, assumed to map from codomain to domain. This method can take out as an optional argument but this is not mandatory.

Examples

# example 1: a matrix encapsulated as a (fake) matrix-free linear operator A_SM = StencilMatrix(V, W) AT_SM = A_SM.transpose() A = MatrixFreeLinearOperator(domain=V, codomain=W, dot=lambda v: A_SM @ v, dot_transpose=lambda v: AT_SM @ v)

# example 2: a truly matrix-free linear operator A = MatrixFreeLinearOperator(domain=V, codomain=V, dot=lambda v: 2*v, dot_transpose=lambda v: 2*v)

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property dtype#
dot(v, out=None)[source]#

Apply linear operator to Vector v. Result is written to Vector out, if provided.

toarray()[source]#

Convert to Numpy 2D array.

tosparse()[source]#
transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.