Data Types | |
type | epetra_multivector |
Public Member Functions | |
Constructor Functions | |
type(epetra_multivector) function, public | epetra_multivector (BlockMap, Num_Vectors, zero) |
Epetra_MultiVector constructor conformal to a BlockMap, optionally zero the newly created vector. More... | |
type(epetra_multivector) function, public | epetra_multivector (this) |
Epetra_MultiVector copy constructor. More... | |
Private Member Functions | |
Post-construction modification routines | |
subroutine | replaceglobalvalue (this, GlobalRow, VectorIndex, ScalarValue, err) |
Replaces value at location (GlobalRow,VectorIndex) with ScalarValue More... | |
subroutine | replaceglobalvalue (this, GlobalBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err) |
Replaces value at location (GlobalBlockRow,BlockRowOffset,VectorIndex) with ScalarValue More... | |
subroutine | replacemyvalue (this, MyRow, VectorIndex, ScalarValue, err) |
Replaces value at location (MyRow,VectorIndex) with ScalarValue More... | |
subroutine | replacemyvalue (this, MyBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err) |
Replaces value at location (MyBlockRow,BlockRowOffset,VectorIndex) with ScalarValue More... | |
subroutine | sumintoglobalvalue (this, GlobalRow, VectorIndex, ScalarValue, err) |
Adds ScalarValue to value at location (GlobalRow,VectorIndex) More... | |
subroutine | sumintoglobalvalue (this, GlobalBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err) |
Adds ScalarValue to value at location (GlobalBlockRow,BlockRowOffset,VectorIndex) More... | |
subroutine | sumintomyvalue (this, MyRow, VectorIndex, ScalarValue, err) |
Adds ScalarValue to value at location (MyRow,VectorIndex) More... | |
subroutine | sumintomyvalue (this, MyBlockRow, BlockRowOffset, VectorIndex, ScalarValue, err) |
Adds ScalarValue to value at location (MyBlockRow,BlockRowOffset,VectorIndex) More... | |
subroutine | putscalar (this, scalar, err) |
Replaces all entries with scalar value. More... | |
subroutine | random (this, err) |
Replaces all entries with random values. More... | |
Mathematical methods | |
real(c_double) function, dimension(:), allocatable | dot (this, x, err) |
Computes the scalar product of corresponding pairs of vectors. More... | |
subroutine | abs (this, A, err) |
Replaces target with element-wise absolute value of input More... | |
subroutine | reciprocal (this, A, err) |
Reciprocal replaces target with element-wise reciprocal value of input More... | |
subroutine | scale (this, scalar_value, err) |
Scales current values this = scalar_value*this More... | |
subroutine | scale (this, scalar_value, MultiVector, err) |
Replaces current values with scaled input this = scalar_value*MultiVector More... | |
real(c_double) function, dimension(:), allocatable | norm1 (this, err) |
Computes 1-norm of each vector in the input More... | |
real(c_double) function, dimension(:), allocatable | norm2 (this, err) |
Computes 2-norm of each vector in the input More... | |
real(c_double) function, dimension(:), allocatable | norminf (this, err) |
Computes infinity norm of each vector in the input More... | |
real(c_double) function, dimension(:), allocatable | normweighted (this, weights, err) |
Computes weighted norm (RMS norm) of each vector in the input More... | |
real(c_double) function, dimension(:), allocatable | minvalue (this, err) |
Computes minimum value of each vector in the input More... | |
real(c_double) function, dimension(:), allocatable | maxvalue (this, err) |
MaxValue: compute maximum value of each vector in the input More... | |
real(c_double) function, dimension(:), allocatable | meanvalue (this, err) |
Computes mean (average) value of each vector in the input More... | |
subroutine | multiply (this, TransA, TransB, ScalarAB, A, B, ScalarThis, err) |
Matrix-matrix multiplication this = ScalarThis*This + ScalarAB*MATMUL(A,B) More... | |
subroutine | multiply (this, ScalarAB, A, B, ScalarThis, err) |
Element-by-element multiplication this = ScalarThis*This + ScalarAB*A*B More... | |
subroutine | reciprocalmultiply (this, ScalarAB, A, B, ScalarThis, err) |
Element-by-element multiplication by reciprocal this = ScalarThis*This + ScalarAB*A/B More... | |
subroutine | update (this, scalarA, A, scalarThis, err) |
Updates with scaled copy of input this = ScalarThis*This + ScalarA*A More... | |
subroutine | update (this, scalarA, A, scalarB, B, scalarThis, err) |
Updates with scaled copies of inputs this = ScalarThis*This + ScalarA*A + ScalarB*B More... | |
Extraction methods | |
real(c_double) function, dimension(:,:), allocatable | extractcopy (this, MyLDA, err) |
Copies multivector contents into target More... | |
Attribute access | |
integer(c_int) function | numvectors (this) |
Number of vectors in multivector More... | |
integer(c_int) function | mylength (this) |
Local vector length More... | |
integer(c_int) function | globallength (this) |
Global vector length More... | |
integer(c_int) function | stride (this) |
Stride between successive vectors in multivector (only meaningful if ConstantStride()==true) More... | |
integer(ft_boolean_t) function | constantstride (this) |
True if stride between successive vectors is constant More... | |
|
private |
Replaces target with element-wise absolute value of input
[in] | a | In The source multivector. |
[out] | err | Returns error information. |
|
private |
True if stride between successive vectors is constant
|
private |
Computes the scalar product of corresponding pairs of vectors.
[in] | x | In The multivector to be used in conjunction with the "this" multivector |
[out] | err | Returns error information. |
type(epetra_multivector) function, public fepetra_multivector::epetra_multivector | ( | class(epetra_multivector), intent(in) | this | ) |
Epetra_MultiVector copy constructor.
type(epetra_multivector) function, public fepetra_multivector::epetra_multivector | ( | class(epetra_blockmap), intent(in) | BlockMap, |
integer(c_int), intent(in) | Num_Vectors, | ||
logical, intent(in) | zero | ||
) |
Epetra_MultiVector constructor conformal to a BlockMap, optionally zero the newly created vector.
[in] | blockmap | In The map to which the vector will conform |
[in] | num_vectors | In Number of vectors in multivector |
[in] | zero | In Optionally zero out the output. |
|
private |
Copies multivector contents into target
The input argument MyLDA is a user request for the size of the output; the actual size will be the maximum between this and the stride of the multivector object.
[in] | mylda | In: Minimum leading dimension of result. |
[out] | err | Returns error information. |
|
private |
Global vector length
|
private |
MaxValue: compute maximum value of each vector in the input
[out] | err | Returns error information. |
|
private |
Computes mean (average) value of each vector in the input
[out] | err | Returns error information. |
|
private |
Computes minimum value of each vector in the input
[out] | err | Returns error information. |
|
private |
Matrix-matrix multiplication this = ScalarThis*This + ScalarAB*MATMUL(A,B)
This routine will compute the product of the two multivectors A and B and use it to update "this", like the BLAS routine GEMM.
[in] | transa | In: Choose transpose status of A |
[in] | transb | In: Choose transpose status of B |
[in] | a | In: Input multivector A |
[in] | b | In: Input multivector B |
[in] | scalarab | In: scale factor for product MATMUL(A,B) |
[in] | scalarthis | In: scale factor for "this" |
[out] | err | Returns error information. |
|
private |
Element-by-element multiplication this = ScalarThis*This + ScalarAB*A*B
[in] | a | In: Input multivector A |
[in] | b | In: Input multivector B |
[in] | scalarab | In: scale factor for product A*B |
[in] | scalarthis | In: scale factor for "this" |
[out] | err | Returns error information. |
|
private |
Local vector length
|
private |
Computes 1-norm of each vector in the input
[out] | err | Returns error information. |
|
private |
Computes 2-norm of each vector in the input
[out] | err | Returns error information. |
|
private |
Computes infinity norm of each vector in the input
[out] | err | Returns error information. |
|
private |
Computes weighted norm (RMS norm) of each vector in the input
[in] | weights | In: the weights. |
[out] | err | Returns error information. |
|
private |
Number of vectors in multivector
|
private |
Replaces all entries with scalar value.
[in] | scalar | In The scalar to which all entries will be set |
[out] | err | Returns error information. |
|
private |
Replaces all entries with random values.
[out] | err | Returns error information. |
|
private |
Reciprocal replaces target with element-wise reciprocal value of input
[in] | a | In The source multivector. |
[out] | err | Returns error information. |
|
private |
Element-by-element multiplication by reciprocal this = ScalarThis*This + ScalarAB*A/B
[in] | a | In: Input multivector A |
[in] | b | In: Input multivector B |
[in] | scalarab | In: scale factor for reciprocal product A/B |
[in] | scalarthis | In: scale factor for "this" |
[out] | err | Returns error information. |
|
private |
Replaces value at location (GlobalBlockRow,BlockRowOffset,VectorIndex) with ScalarValue
The index of the specified location must correspond to an index owned by the map on the calling process; i.e. no communication takes place.
[in] | globalblockrow | In The global block row to be set |
[in] | blockrowoffset | In The global block row offest to be set |
[in] | vectorindex | In The vector index within the multivector to be set |
[in] | scalarvalue | In The scalar value to be set |
[out] | err | Returns error information. |
|
private |
Replaces value at location (GlobalRow,VectorIndex) with ScalarValue
The index of the specified location must correspond to an index owned by the map on the calling process; i.e. no communication takes place.
[in] | globalrow | In The global row to be set |
[in] | vectorindex | In The vector index within the multivector to be set |
[in] | scalarvalue | In The value to be set |
[out] | err | Returns error information. |
|
private |
Replaces value at location (MyBlockRow,BlockRowOffset,VectorIndex) with ScalarValue
[in] | myblockrow | In The local block row to be set |
[in] | blockrowoffset | In The local block row offset to be set |
[in] | vectorindex | In The vector index within the multivector to be set |
[in] | scalarvalue | In The scalar value to be set |
|
private |
Replaces value at location (MyRow,VectorIndex) with ScalarValue
The index of the specified location must be that of a locally owned element.
[in] | myrow | In The local row to be set |
[in] | vectorindex | In The vector index within the multivector to be set |
[in] | scalarvalue | In The scalar value to be set |
[out] | err | Returns error information. |
|
private |
Scales current values this = scalar_value*this
[in] | scalar_value | In The scale factor. |
[out] | err | Returns error information. |
|
private |
Replaces current values with scaled input this = scalar_value*MultiVector
[in] | multivector | In The source multivector. |
[in] | scalar_value | In The scale factor. |
[out] | err | Returns error information. |
|
private |
Stride between successive vectors in multivector (only meaningful if ConstantStride()==true)
|
private |
Adds ScalarValue to value at location (GlobalRow,VectorIndex)
[in] | globalrow | In The global row to be modified |
[in] | vectorindex | In The vector index within the multivector to be modified |
[in] | scalarvalue | In The scalar value to be added |
|
private |
Adds ScalarValue to value at location (GlobalBlockRow,BlockRowOffset,VectorIndex)
[in] | globalblockrow | In The global block row to be modified |
[in] | blockrowoffset | In The global block row offset to be modified |
[in] | vectorindex | In The vector index within the multivector to be modified |
[in] | scalarvalue | In The scalar value to be added |
[out] | err | Returns error information. |
|
private |
Adds ScalarValue to value at location (MyBlockRow,BlockRowOffset,VectorIndex)
[in] | myblockrow | In The local block row to be modified |
[in] | blockrowoffset | In The local block row offset to be modified |
[in] | vectorindex | In The vector index within the multivector to be modified |
[in] | scalarvalue | In The scalar value to be added |
[out] | err | Returns error information. |
|
private |
Adds ScalarValue to value at location (MyRow,VectorIndex)
[in] | myrow | In The local row to be modified |
[in] | vectorindex | In The vector index within the multivector to be modified |
[in] | scalarvalue | In The scalar value to be added |
[out] | err | Returns error information. |
|
private |
Updates with scaled copy of input this = ScalarThis*This + ScalarA*A
[in] | a | In: Input multivector A |
[in] | scalara | In: scale factor for p A |
[in] | scalarthis | In: scale factor for "this" |
[out] | err | Returns error information. |
|
private |
Updates with scaled copies of inputs this = ScalarThis*This + ScalarA*A + ScalarB*B
[in] | scalara | In: scale factor for A |
[in] | a | In: Input multivector A |
[in] | scalarb | In: scale factor for B |
[in] | b | In: Input multivector B |
[in] | scalarthis | In: scale factor for "this" |
[out] | err | Returns error information. |