pymor.vectorarrays.numpy¶
Module Contents¶
Classes¶
|
|
|
- class pymor.vectorarrays.numpy.NumpyVectorArray(array, space)[source]¶
Bases:
pymor.vectorarrays.interface.VectorArrayVectorArrayimplementation viaNumPy arrays.This is the default
VectorArraytype used by allOperatorsin pyMOR’s discretization toolkit. Moreover, all reducedOperatorsare based onNumpyVectorArray.This class is just a thin wrapper around the underlying
NumPy array. Thus, while operations likeaxpyorinnerwill be quite efficient, removing or appending vectors will be costly.Warning
This class is not intended to be instantiated directly. Use the associated
VectorSpaceinstead.- 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 returnedNumPy arraymight alter the originalVectorArray. IfTruealways a copy of the array data is made.
- __getitem__(self, ind)[source]¶
Return a
VectorArrayview onto a subset of the vectors in the array.
- copy(self, deep=False, *, _ind=None)[source]¶
Returns a copy of the array.
All
VectorArrayimplementations in pyMOR have copy-on-write semantics: if not specified otherwise by settingdeeptoTrue, 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
VectorArraycontaining the vectors to be appended.- remove_from_other
If
True, the appended vectors are removed fromother. 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
alphais a scalar, each vector is multiplied by this scalar. Otherwise,alphahas to be a one-dimensionalNumPy arrayof the same length asselfcontaining the factors for each vector.Parameters
- alpha
The scalar coefficient or one-dimensional
NumPy arrayof coefficients with which the vectors inselfare multiplied.
- axpy(self, alpha, x, *, _ind=None)[source]¶
BLAS AXPY operation.
This method forms the sum
self = alpha*x + self
If the length of
xis 1, the samexvector is used for all vectors inself. Otherwise, the lengths ofselfandxhave to agree. Ifalphais a scalar, eachxvector is multiplied with the same factoralpha. Otherwise,alphahas to be a one-dimensionalNumPy arrayof the same length asselfcontaining the coefficients for eachxvector.Parameters
- alpha
The scalar coefficient or one-dimensional
NumPy arrayof coefficients with which the vectors inxare multiplied.- x
A
VectorArraycontaining the x-summands.
- inner(self, other, product=None, *, _ind=None)[source]¶
Returns the inner products between
VectorArrayelements.If
productisNone, the Euclidean inner product between thedofsofselfandotherare returned, i.e.U.inner(V)
is equivalent to:
U.dofs(np.arange(U.dim)) @ V.dofs(np.arange(V.dim)).T
(Note, that
dofsis only intended to be called for a small number of DOF indices.)If a
productOperatoris specified, thisOperatoris used to compute the inner products usingapply2, 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
VectorArraycontaining the second factors.- product
If not
NoneanOperatorrepresenting the inner product bilinear form.
- pairwise_inner(self, other, product=None, *, _ind=None)[source]¶
Returns the pairwise inner products between
VectorArrayelements.If
productisNone, the Euclidean inner product between thedofsofselfandotherare 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
dofsis only intended to be called for a small number of DOF indices.)If a
productOperatoris specified, thisOperatoris used to compute the inner products usingpairwise_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
VectorArraycontaining the second factors.- product
If not
NoneanOperatorrepresenting the inner product bilinear form.
- lincomb(self, coefficients, *, _ind=None)[source]¶
Returns linear combinations of the vectors contained in the array.
Parameters
- coefficients
A
NumPy arrayof dimension 1 or 2 containing the linear coefficients.coefficients.shape[-1]has to agree withlen(self).
Returns
A
VectorArrayresultsuch thatresult[i] = ∑ self[j] * coefficients[i,j]
in case
coefficientsis of dimension 2, otherwiselen(result) == 1andresult[0] = ∑ self[j] * coefficients[j].
- sup_norm(self, *, _ind=None)[source]¶
The l-infinity-norms of the vectors contained in the array.
Returns
A
NumPy arrayresultsuch thatresult[i]contains the norm ofself[i].
- dofs(self, dof_indices, *, _ind=None)[source]¶
Extract DOFs of the vectors contained in the array.
Parameters
- dof_indices
List or 1D
NumPy arrayof indices of the DOFs that are to be returned.
Returns
A
NumPy arrayresultsuch thatresult[i, j]is thedof_indices[j]-th DOF of thei-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 arraycontaining for each vector a DOF index at which the maximum is attained.- max_val
NumPy arraycontaining for each vector the maximum absolute value of its DOFs.
- class pymor.vectorarrays.numpy.NumpyVectorSpace(dim, id=None)[source]¶
Bases:
pymor.vectorarrays.interface.VectorSpaceVectorSpaceofNumpyVectorArrays.Parameters
- dim
The dimension of the vectors contained in the space.
- id
See
id.
- zeros(self, count=1, reserve=0)[source]¶
Create a
VectorArrayof null vectorsParameters
- count
The number of vectors.
- reserve
Hint for the backend to which length the array will grow.
Returns
A
VectorArraycontainingcountvectors with each component zero.
- full(self, value, count=1, reserve=0)[source]¶
Create a
VectorArrayof 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
VectorArraycontainingcountvectors with each DOF set tovalue.
- random(self, count=1, distribution='uniform', random_state=None, seed=None, reserve=0, **kwargs)[source]¶
Create a
VectorArrayof 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
VectorSpaceimplementations.Parameters
- count
The number of vectors.
- distribution
Random distribution to use (
'uniform','normal').- low
Lower bound for
'uniform'distribution (defaults to0).- high
Upper bound for
'uniform'distribution (defaults to1).- loc
Mean for
'normal'distribution (defaults to0).- scale
Standard deviation for
'normal'distribution (defaults to1).- random_state
RandomStateto use for sampling. IfNone, a new random state is generated usingseedas random seed, or thedefaultrandom 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
VectorArrayfrom raw data.This method is used in the implementation of
OperatorsandModelsto create newVectorArraysfrom 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
VectorArrayfrom aNumPy arrayNote that this method will not be supported by all vector space implementations.
Parameters
- data
NumPyarray of shape(len, dim)wherelenis the number of vectors anddimtheir dimension.- ensure_copy
If
False, modifying the returnedVectorArraymight alter the originalNumPy array. IfTruealways a copy of the array data is made.
Returns
A
VectorArraywithdataas data.
- class pymor.vectorarrays.numpy.NumpyVectorArrayView(array, ind)[source]¶
Bases:
NumpyVectorArrayVectorArrayimplementation viaNumPy arrays.This is the default
VectorArraytype used by allOperatorsin pyMOR’s discretization toolkit. Moreover, all reducedOperatorsare based onNumpyVectorArray.This class is just a thin wrapper around the underlying
NumPy array. Thus, while operations likeaxpyorinnerwill be quite efficient, removing or appending vectors will be costly.Warning
This class is not intended to be instantiated directly. Use the associated
VectorSpaceinstead.- 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 returnedNumPy arraymight alter the originalVectorArray. IfTruealways a copy of the array data is made.
- __getitem__(self, ind)[source]¶
Return a
VectorArrayview onto a subset of the vectors in the array.
- append(self, other, remove_from_other=False)[source]¶
Append vectors to the array.
Parameters
- other
A
VectorArraycontaining the vectors to be appended.- remove_from_other
If
True, the appended vectors are removed fromother. 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
VectorArrayimplementations in pyMOR have copy-on-write semantics: if not specified otherwise by settingdeeptoTrue, 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
alphais a scalar, each vector is multiplied by this scalar. Otherwise,alphahas to be a one-dimensionalNumPy arrayof the same length asselfcontaining the factors for each vector.Parameters
- alpha
The scalar coefficient or one-dimensional
NumPy arrayof coefficients with which the vectors inselfare multiplied.
- axpy(self, alpha, x)[source]¶
BLAS AXPY operation.
This method forms the sum
self = alpha*x + self
If the length of
xis 1, the samexvector is used for all vectors inself. Otherwise, the lengths ofselfandxhave to agree. Ifalphais a scalar, eachxvector is multiplied with the same factoralpha. Otherwise,alphahas to be a one-dimensionalNumPy arrayof the same length asselfcontaining the coefficients for eachxvector.Parameters
- alpha
The scalar coefficient or one-dimensional
NumPy arrayof coefficients with which the vectors inxare multiplied.- x
A
VectorArraycontaining the x-summands.
- inner(self, other, product=None)[source]¶
Returns the inner products between
VectorArrayelements.If
productisNone, the Euclidean inner product between thedofsofselfandotherare returned, i.e.U.inner(V)
is equivalent to:
U.dofs(np.arange(U.dim)) @ V.dofs(np.arange(V.dim)).T
(Note, that
dofsis only intended to be called for a small number of DOF indices.)If a
productOperatoris specified, thisOperatoris used to compute the inner products usingapply2, 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
VectorArraycontaining the second factors.- product
If not
NoneanOperatorrepresenting the inner product bilinear form.
- pairwise_inner(self, other, product=None)[source]¶
Returns the pairwise inner products between
VectorArrayelements.If
productisNone, the Euclidean inner product between thedofsofselfandotherare 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
dofsis only intended to be called for a small number of DOF indices.)If a
productOperatoris specified, thisOperatoris used to compute the inner products usingpairwise_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
VectorArraycontaining the second factors.- product
If not
NoneanOperatorrepresenting the inner product bilinear form.
- lincomb(self, coefficients)[source]¶
Returns linear combinations of the vectors contained in the array.
Parameters
- coefficients
A
NumPy arrayof dimension 1 or 2 containing the linear coefficients.coefficients.shape[-1]has to agree withlen(self).
Returns
A
VectorArrayresultsuch thatresult[i] = ∑ self[j] * coefficients[i,j]
in case
coefficientsis of dimension 2, otherwiselen(result) == 1andresult[0] = ∑ self[j] * coefficients[j].
- sup_norm(self)[source]¶
The l-infinity-norms of the vectors contained in the array.
Returns
A
NumPy arrayresultsuch thatresult[i]contains the norm ofself[i].
- dofs(self, dof_indices)[source]¶
Extract DOFs of the vectors contained in the array.
Parameters
- dof_indices
List or 1D
NumPy arrayof indices of the DOFs that are to be returned.
Returns
A
NumPy arrayresultsuch thatresult[i, j]is thedof_indices[j]-th DOF of thei-th vector of the array.
- amax(self)[source]¶
The maximum absolute value of the DOFs contained in the array.
Returns
- max_ind
NumPy arraycontaining for each vector a DOF index at which the maximum is attained.- max_val
NumPy arraycontaining for each vector the maximum absolute value of its DOFs.