pymor.vectorarrays.numpy

Module Contents

Classes

NumpyVectorArray

VectorArray implementation via NumPy arrays.

NumpyVectorSpace

VectorSpace of NumpyVectorArrays.

NumpyVectorArrayView

VectorArray implementation via NumPy arrays.

class pymor.vectorarrays.numpy.NumpyVectorArray(array, space)[source]

Bases: pymor.vectorarrays.interface.VectorArray

VectorArray implementation via NumPy arrays.

This is the default VectorArray type used by all Operators in pyMOR’s discretization toolkit. Moreover, all reduced Operators are based on NumpyVectorArray.

This class is just a thin wrapper around the underlying NumPy array. Thus, while operations like axpy or inner will be quite efficient, removing or appending vectors will be costly.

Warning

This class is not intended to be instantiated directly. Use the associated VectorSpace instead.

__radd__[source]
to_numpy(self, ensure_copy=False)[source]

Return (len(self), self.dim) NumPy Array with the data stored in the array.

Parameters

ensure_copy

If False, modifying the returned NumPy array might alter the original VectorArray. If True always a copy of the array data is made.

property _data(self)[source]

Return NumPy Array view on data for hacking / interactive use.

property real(self)[source]

Real part.

property imag(self)[source]

Imaginary part.

conj(self)[source]

Complex conjugation.

__len__(self)[source]

The number of vectors in the array.

__getitem__(self, ind)[source]

Return a VectorArray view onto a subset of the vectors in the array.

__delitem__(self, ind)[source]

Remove vectors from the array.

copy(self, deep=False, *, _ind=None)[source]

Returns a copy of the array.

All VectorArray implementations in pyMOR have copy-on-write semantics: if not specified otherwise by setting deep to True, the returned copy will hold a handle to the same array data as the original array, and a deep copy of the data will only be performed when one of the arrays is modified.

Note that for NumpyVectorArray, a deep copy is always performed when only some vectors in the array are copied.

Parameters

deep

Ensure that an actual copy of the array data is made (see above).

Returns

A copy of the VectorArray.

append(self, other, remove_from_other=False)[source]

Append vectors to the array.

Parameters

other

A VectorArray containing the vectors to be appended.

remove_from_other

If True, the appended vectors are removed from other. For list-like implementations this can be used to prevent unnecessary copies of the involved vectors.

scal(self, alpha, *, _ind=None)[source]

BLAS SCAL operation (in-place scalar multiplication).

This method calculates

self = alpha*self

If alpha is a scalar, each vector is multiplied by this scalar. Otherwise, alpha has to be a one-dimensional NumPy array of the same length as self containing the factors for each vector.

Parameters

alpha

The scalar coefficient or one-dimensional NumPy array of coefficients with which the vectors in self are multiplied.

axpy(self, alpha, x, *, _ind=None)[source]

BLAS AXPY operation.

This method forms the sum

self = alpha*x + self

If the length of x is 1, the same x vector is used for all vectors in self. Otherwise, the lengths of self and x have to agree. If alpha is a scalar, each x vector is multiplied with the same factor alpha. Otherwise, alpha has to be a one-dimensional NumPy array of the same length as self containing the coefficients for each x vector.

Parameters

alpha

The scalar coefficient or one-dimensional NumPy array of coefficients with which the vectors in x are multiplied.

x

A VectorArray containing the x-summands.

inner(self, other, product=None, *, _ind=None)[source]

Returns the inner products between VectorArray elements.

If product is None, the Euclidean inner product between the dofs of self and other are returned, i.e.

U.inner(V)

is equivalent to:

U.dofs(np.arange(U.dim)) @ V.dofs(np.arange(V.dim)).T

(Note, that dofs is only intended to be called for a small number of DOF indices.)

If a product Operator is specified, this Operator is used to compute the inner products using apply2, i.e. U.inner(V, product) is equivalent to:

product.apply2(U, V)

which in turn is, by default, implemented as:

U.inner(product.apply(V))

In the case of complex numbers, this is antilinear in the first argument, i.e. in ‘self’. Complex conjugation is done in the first argument because most numerical software in the community handles it this way: Numpy, DUNE, FEniCS, Eigen, Matlab and BLAS do complex conjugation in the first argument, only PetSc and deal.ii do complex conjugation in the second argument.

Parameters

other

A VectorArray containing the second factors.

product

If not None an Operator representing the inner product bilinear form.

Returns

A NumPy array result such that

result[i, j] = ( self[i], other[j] ).

pairwise_inner(self, other, product=None, *, _ind=None)[source]

Returns the pairwise inner products between VectorArray elements.

If product is None, the Euclidean inner product between the dofs of self and other are returned, i.e.

U.pairwise_inner(V)

is equivalent to:

np.sum(U.dofs(np.arange(U.dim)) * V.dofs(np.arange(V.dim)), axis=-1)

(Note, that dofs is only intended to be called for a small number of DOF indices.)

If a product Operator is specified, this Operator is used to compute the inner products using pairwise_apply2, i.e. U.inner(V, product) is equivalent to:

product.pairwise_apply2(U, V)

which in turn is, by default, implemented as:

U.pairwise_inner(product.apply(V))

In the case of complex numbers, this is antilinear in the first argument, i.e. in ‘self’. Complex conjugation is done in the first argument because most numerical software in the community handles it this way: Numpy, DUNE, FEniCS, Eigen, Matlab and BLAS do complex conjugation in the first argument, only PetSc and deal.ii do complex conjugation in the second argument.

Parameters

other

A VectorArray containing the second factors.

product

If not None an Operator representing the inner product bilinear form.

Returns

A NumPy array result such that

result[i] = ( self[i], other[i] ).

lincomb(self, coefficients, *, _ind=None)[source]

Returns linear combinations of the vectors contained in the array.

Parameters

coefficients

A NumPy array of dimension 1 or 2 containing the linear coefficients. coefficients.shape[-1] has to agree with len(self).

Returns

A VectorArray result such that

result[i] = ∑ self[j] * coefficients[i,j]

in case coefficients is of dimension 2, otherwise len(result) == 1 and

result[0] = ∑ self[j] * coefficients[j].

_norm(self, *, _ind=None)[source]

Implementation of norm for the case that no product is given.

_norm2(self, *, _ind=None)[source]

Implementation of norm2 for the case that no product is given.

sup_norm(self, *, _ind=None)[source]

The l-infinity-norms of the vectors contained in the array.

Returns

A NumPy array result such that result[i] contains the norm of self[i].

dofs(self, dof_indices, *, _ind=None)[source]

Extract DOFs of the vectors contained in the array.

Parameters

dof_indices

List or 1D NumPy array of indices of the DOFs that are to be returned.

Returns

A NumPy array result such that result[i, j] is the dof_indices[j]-th DOF of the i-th vector of the array.

amax(self, *, _ind=None)[source]

The maximum absolute value of the DOFs contained in the array.

Returns

max_ind

NumPy array containing for each vector a DOF index at which the maximum is attained.

max_val

NumPy array containing for each vector the maximum absolute value of its DOFs.

__str__(self)[source]

Return str(self).

_format_repr(self, max_width, verbosity)[source]
__del__(self)[source]
_deep_copy(self)[source]
__add__(self, other)[source]
__iadd__(self, other)[source]
__sub__(self, other)[source]
__isub__(self, other)[source]
__mul__(self, other)[source]
__imul__(self, other)[source]
__neg__(self)[source]
class pymor.vectorarrays.numpy.NumpyVectorSpace(dim, id=None)[source]

Bases: pymor.vectorarrays.interface.VectorSpace

VectorSpace of NumpyVectorArrays.

Parameters

dim

The dimension of the vectors contained in the space.

id

See id.

__eq__(self, other)[source]

Return self==value.

__hash__(self)[source]

Return hash(self).

zeros(self, count=1, reserve=0)[source]

Create a VectorArray of null vectors

Parameters

count

The number of vectors.

reserve

Hint for the backend to which length the array will grow.

Returns

A VectorArray containing count vectors with each component zero.

full(self, value, count=1, reserve=0)[source]

Create a VectorArray of vectors with all DOFs set to the same value.

Parameters

value

The value each DOF should be set to.

count

The number of vectors.

reserve

Hint for the backend to which length the array will grow.

Returns

A VectorArray containing count vectors with each DOF set to value.

random(self, count=1, distribution='uniform', random_state=None, seed=None, reserve=0, **kwargs)[source]

Create a VectorArray of vectors with random entries.

Supported random distributions:

'uniform': Uniform distribution in half-open interval
           [`low`, `high`).
'normal':  Normal (Gaussian) distribution with mean
           `loc` and standard deviation `scale`.

Note that not all random distributions are necessarily implemented by all VectorSpace implementations.

Parameters

count

The number of vectors.

distribution

Random distribution to use ('uniform', 'normal').

low

Lower bound for 'uniform' distribution (defaults to 0).

high

Upper bound for 'uniform' distribution (defaults to 1).

loc

Mean for 'normal' distribution (defaults to 0).

scale

Standard deviation for 'normal' distribution (defaults to 1).

random_state

RandomState to use for sampling. If None, a new random state is generated using seed as random seed, or the default random state is used.

seed

If not None, a new random state with this seed is used.

reserve

Hint for the backend to which length the array will grow.

make_array(cls, obj, id=None)[source]

Create a VectorArray from raw data.

This method is used in the implementation of Operators and Models to create new VectorArrays from raw data of the underlying solver backends. The ownership of the data is transferred to the newly created array.

The exact signature of this method depends on the wrapped solver backend.

from_numpy(cls, data, id=None, ensure_copy=False)[source]

Create a VectorArray from a NumPy array

Note that this method will not be supported by all vector space implementations.

Parameters

data

NumPy array of shape (len, dim) where len is the number of vectors and dim their dimension.

ensure_copy

If False, modifying the returned VectorArray might alter the original NumPy array. If True always a copy of the array data is made.

Returns

A VectorArray with data as data.

from_file(cls, path, key=None, single_vector=False, transpose=False, id=None)[source]
classmethod _array_factory(cls, array, space=None, id=None)[source]
property is_scalar(self)[source]
class pymor.vectorarrays.numpy.NumpyVectorArrayView(array, ind)[source]

Bases: NumpyVectorArray

VectorArray implementation via NumPy arrays.

This is the default VectorArray type used by all Operators in pyMOR’s discretization toolkit. Moreover, all reduced Operators are based on NumpyVectorArray.

This class is just a thin wrapper around the underlying NumPy array. Thus, while operations like axpy or inner will be quite efficient, removing or appending vectors will be costly.

Warning

This class is not intended to be instantiated directly. Use the associated VectorSpace instead.

is_view = True[source]
__radd__[source]
to_numpy(self, ensure_copy=False)[source]

Return (len(self), self.dim) NumPy Array with the data stored in the array.

Parameters

ensure_copy

If False, modifying the returned NumPy array might alter the original VectorArray. If True always a copy of the array data is made.

__len__(self)[source]

The number of vectors in the array.

__getitem__(self, ind)[source]

Return a VectorArray view onto a subset of the vectors in the array.

__delitem__(self)[source]

Remove vectors from the array.

append(self, other, remove_from_other=False)[source]

Append vectors to the array.

Parameters

other

A VectorArray containing the vectors to be appended.

remove_from_other

If True, the appended vectors are removed from other. For list-like implementations this can be used to prevent unnecessary copies of the involved vectors.

copy(self, deep=False)[source]

Returns a copy of the array.

All VectorArray implementations in pyMOR have copy-on-write semantics: if not specified otherwise by setting deep to True, the returned copy will hold a handle to the same array data as the original array, and a deep copy of the data will only be performed when one of the arrays is modified.

Note that for NumpyVectorArray, a deep copy is always performed when only some vectors in the array are copied.

Parameters

deep

Ensure that an actual copy of the array data is made (see above).

Returns

A copy of the VectorArray.

scal(self, alpha)[source]

BLAS SCAL operation (in-place scalar multiplication).

This method calculates

self = alpha*self

If alpha is a scalar, each vector is multiplied by this scalar. Otherwise, alpha has to be a one-dimensional NumPy array of the same length as self containing the factors for each vector.

Parameters

alpha

The scalar coefficient or one-dimensional NumPy array of coefficients with which the vectors in self are multiplied.

axpy(self, alpha, x)[source]

BLAS AXPY operation.

This method forms the sum

self = alpha*x + self

If the length of x is 1, the same x vector is used for all vectors in self. Otherwise, the lengths of self and x have to agree. If alpha is a scalar, each x vector is multiplied with the same factor alpha. Otherwise, alpha has to be a one-dimensional NumPy array of the same length as self containing the coefficients for each x vector.

Parameters

alpha

The scalar coefficient or one-dimensional NumPy array of coefficients with which the vectors in x are multiplied.

x

A VectorArray containing the x-summands.

inner(self, other, product=None)[source]

Returns the inner products between VectorArray elements.

If product is None, the Euclidean inner product between the dofs of self and other are returned, i.e.

U.inner(V)

is equivalent to:

U.dofs(np.arange(U.dim)) @ V.dofs(np.arange(V.dim)).T

(Note, that dofs is only intended to be called for a small number of DOF indices.)

If a product Operator is specified, this Operator is used to compute the inner products using apply2, i.e. U.inner(V, product) is equivalent to:

product.apply2(U, V)

which in turn is, by default, implemented as:

U.inner(product.apply(V))

In the case of complex numbers, this is antilinear in the first argument, i.e. in ‘self’. Complex conjugation is done in the first argument because most numerical software in the community handles it this way: Numpy, DUNE, FEniCS, Eigen, Matlab and BLAS do complex conjugation in the first argument, only PetSc and deal.ii do complex conjugation in the second argument.

Parameters

other

A VectorArray containing the second factors.

product

If not None an Operator representing the inner product bilinear form.

Returns

A NumPy array result such that

result[i, j] = ( self[i], other[j] ).

pairwise_inner(self, other, product=None)[source]

Returns the pairwise inner products between VectorArray elements.

If product is None, the Euclidean inner product between the dofs of self and other are returned, i.e.

U.pairwise_inner(V)

is equivalent to:

np.sum(U.dofs(np.arange(U.dim)) * V.dofs(np.arange(V.dim)), axis=-1)

(Note, that dofs is only intended to be called for a small number of DOF indices.)

If a product Operator is specified, this Operator is used to compute the inner products using pairwise_apply2, i.e. U.inner(V, product) is equivalent to:

product.pairwise_apply2(U, V)

which in turn is, by default, implemented as:

U.pairwise_inner(product.apply(V))

In the case of complex numbers, this is antilinear in the first argument, i.e. in ‘self’. Complex conjugation is done in the first argument because most numerical software in the community handles it this way: Numpy, DUNE, FEniCS, Eigen, Matlab and BLAS do complex conjugation in the first argument, only PetSc and deal.ii do complex conjugation in the second argument.

Parameters

other

A VectorArray containing the second factors.

product

If not None an Operator representing the inner product bilinear form.

Returns

A NumPy array result such that

result[i] = ( self[i], other[i] ).

lincomb(self, coefficients)[source]

Returns linear combinations of the vectors contained in the array.

Parameters

coefficients

A NumPy array of dimension 1 or 2 containing the linear coefficients. coefficients.shape[-1] has to agree with len(self).

Returns

A VectorArray result such that

result[i] = ∑ self[j] * coefficients[i,j]

in case coefficients is of dimension 2, otherwise len(result) == 1 and

result[0] = ∑ self[j] * coefficients[j].

_norm(self)[source]

Implementation of norm for the case that no product is given.

_norm2(self)[source]

Implementation of norm2 for the case that no product is given.

sup_norm(self)[source]

The l-infinity-norms of the vectors contained in the array.

Returns

A NumPy array result such that result[i] contains the norm of self[i].

dofs(self, dof_indices)[source]

Extract DOFs of the vectors contained in the array.

Parameters

dof_indices

List or 1D NumPy array of indices of the DOFs that are to be returned.

Returns

A NumPy array result such that result[i, j] is the dof_indices[j]-th DOF of the i-th vector of the array.

amax(self)[source]

The maximum absolute value of the DOFs contained in the array.

Returns

max_ind

NumPy array containing for each vector a DOF index at which the maximum is attained.

max_val

NumPy array containing for each vector the maximum absolute value of its DOFs.

__str__(self)[source]

Return str(self).

__add__(self, other)[source]
__iadd__(self, other)[source]
__sub__(self, other)[source]
__isub__(self, other)[source]
__mul__(self, other)[source]
__imul__(self, other)[source]
__neg__(self)[source]
__del__(self)[source]
__repr__(self)[source]

Return repr(self).