Data Types | |
type | epetra_crsmatrix |
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed row matrices. More... | |
Public Member Functions | |
Constructor Functions | |
type(epetra_crsmatrix) function | epetra_crsmatrix (CV, Row_Map, NumEntriesPerRow, StaticProfile) |
Epetra_CrsMatrix constructor with variable number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage. More... | |
type(epetra_crsmatrix) function | epetra_crsmatrix (CV, Row_Map, NumEntriesPerRow, StaticProfile) |
Epetra_CrsMatrix constructor with fixed number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage. More... | |
type(epetra_crsmatrix) function | epetra_crsmatrix (CV, Row_Map, Col_Map, NumEntriesPerRow, StaticProfile) |
Epetra_CrsMatrix constructor with variable number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage. More... | |
type(epetra_crsmatrix) function | epetra_crsmatrix (CV, Row_Map, Col_Map, NumEntriesPerRow, StaticProfile) |
Epetra_CrsMatrix constructor with fixed number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage. More... | |
type(epetra_crsmatrix) function | epetra_crsmatrix (this) |
Epetra_CrsMatrix copy constructor. Creates a copy of a Epetra_CrsMatrix object. More... | |
Insertion/Replace/SumInto methods | |
subroutine | putscalar (this, scalar, err) |
Initialize all values in the matrix with constant value. More... | |
subroutine | scale (this, ScalarConstant, err) |
Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A). More... | |
subroutine | insertglobalvalues (this, GlobalRow, values, indices, err) |
Insert a list of elements in a given global row of the matrix. This method is used to construct a matrix for the first time. It cannot be used if the matrix structure has already been fixed (via a call to FillComplete()). More... | |
subroutine | replaceglobalvalues (this, GlobalRow, values, indices, err) |
Replace specified existing values with this list of entries for a given global row of the matrix. More... | |
Transformation methods | |
subroutine | fillcomplete (this, OptimizeDataStorage, err) |
Signal that data entry is complete. Perform transformations to local index space. This version of FillComplete assumes that the domain and range distributions are identical to the matrix row distributions. More... | |
subroutine | fillcomplete (this, DomainMap, RangeMap, OptimizeDataStorage, err) |
Signal that data entry is complete. Perform transformations to local index space. This version of FillComplete requires the explicit specification of the domain and range distribution maps. These maps are used for importing and exporting vector and multi-vector elements that are needed for distributed matrix computations. More... | |
Extraction Methods | |
subroutine | extractglobalrowcopy (this, GlobalRow, values, indices, err) |
Returns a copy of the specified global row in user-provided arrays. More... | |
integer(c_int) function | nummyentries (this, MyRow) |
Returns the current number of nonzero entries in specified local row on this processor. More... | |
integer(c_int) function | maxnumentries (this) |
Returns the maximum number of nonzero entries across all rows on this processor. More... | |
Additional methods required to implement Epetra_RowMatrix interface | |
integer(c_int) function | nummyrowentries (this, MyRow) |
Return the current number of values stored for the specified local row. More... | |
Computational Methods | |
subroutine | multiply_vector (this, TransA, x, y, err) |
Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_Vector x in y. More... | |
subroutine | multiply_multivector (this, TransA, x, y, err) |
Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_MultiVector X in Y. More... | |
Local/Global ID method | |
logical function | myglobalrow (this, GID) |
Returns true of GID is owned by the calling processor, otherwise it returns false. More... | |
Attribute Accessor Methods | |
type(epetra_map) function | rowmatrixrowmap (this) |
Returns the Epetra_Map object associated with the rows of this matrix. More... | |
type(epetra_map) function | rowmap (this) |
Returns the Epetra_Map object associated with the rows of this matrix. More... | |
integer(c_int) function | numglobalentries (this, GlobalRow) |
Returns the current number of nonzero entries in specified global row on this processor. More... | |
type(epetra_crsmatrix) function fepetra_crsmatrix::epetra_crsmatrix | ( | integer(ft_epetra_dataaccess_e_t), intent(in) | CV, |
class(epetra_map), intent(in) | Row_Map, | ||
integer(c_int), dimension(:), intent(in) | NumEntriesPerRow, | ||
logical, optional | StaticProfile | ||
) |
Epetra_CrsMatrix constructor with variable number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage.
[in] | cv | (In) An FT_Epetra_DataAccess_E_t enumerated type set to Copy or View. |
[in] | row_map | (In) An Epetra_Map defining the numbering and distribution of matrix rows. |
[in] | numentriesperrow | (In) An integer array of length NumRows such that NumEntriesPerRow(i) indicates the (approximate if StaticProfile=false) number of entries in the ith row. |
staticprofile | (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
type(epetra_crsmatrix) function fepetra_crsmatrix::epetra_crsmatrix | ( | integer(ft_epetra_dataaccess_e_t), intent(in) | CV, |
class(epetra_map), intent(in) | Row_Map, | ||
integer(c_int), intent(in) | NumEntriesPerRow, | ||
logical, optional | StaticProfile | ||
) |
Epetra_CrsMatrix constructor with fixed number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage.
[in] | cv | (In) An FT_Epetra_DataAccess_E_t enumerated type set to Copy or View. |
[in] | row_map | (In) An Epetra_Map defining the numbering and distribution of matrix rows. |
[in] | numentriesperrow | (In) An integer that indicates the (approximate) number of entries in the each row. Note that it is possible to use 0 for this value and let fill occur during the insertion phase. |
staticprofile | (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
type(epetra_crsmatrix) function fepetra_crsmatrix::epetra_crsmatrix | ( | integer(ft_epetra_dataaccess_e_t), intent(in) | CV, |
class(epetra_map), intent(in) | Row_Map, | ||
class(epetra_map), intent(in) | Col_Map, | ||
integer(c_int), intent(in) | NumEntriesPerRow, | ||
logical, optional | StaticProfile | ||
) |
Epetra_CrsMatrix constructor with fixed number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage.
[in] | cv | (In) An FT_Epetra_DataAccess_E_t enumerated type set to Copy or View. |
[in] | row_map | (In) An Epetra_Map defining the numbering and distribution of matrix rows. |
[in] | col_map | (In) An Epetra_Map defining the set of column-indices that appear in each processor's locally owned matrix rows. |
[in] | numentriesperrow | (In) An integer that indicates the (approximate if StaticProfile=false) number of entries in the each row. Note that it is possible to use 0 for this value and let fill occur during the insertion phase. |
staticprofile | (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
type(epetra_crsmatrix) function fepetra_crsmatrix::epetra_crsmatrix | ( | type(epetra_crsmatrix), intent(in) | this | ) |
Epetra_CrsMatrix copy constructor. Creates a copy of a Epetra_CrsMatrix object.
type(epetra_crsmatrix) function fepetra_crsmatrix::epetra_crsmatrix | ( | integer(ft_epetra_dataaccess_e_t), intent(in) | CV, |
class(epetra_map), intent(in) | Row_Map, | ||
class(epetra_map), intent(in) | Col_Map, | ||
integer(c_int), dimension(:), intent(in) | NumEntriesPerRow, | ||
logical, optional | StaticProfile | ||
) |
Epetra_CrsMatrix constructor with variable number of indices per row. Creates a Epetra_CrsMatrix object and allocates storage.
[in] | cv | (In) An FT_Epetra_DataAccess_E_t enumerated type set to Copy or View. |
[in] | row_map | (In) An Epetra_Map defining the numbering and distribution of matrix rows. |
[in] | col_map | (In) An Epetra_Map defining the set of column-indices that appear in each processor's locally owned matrix rows. |
[in] | numentriesperrow | (In) An integer array of length NumRows such that NumEntriesPerRow(i) indicates the (approximate if StaticProfile=false) number of entries in the ith row. |
staticprofile | (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
subroutine fepetra_crsmatrix::extractglobalrowcopy | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | GlobalRow, | ||
real(c_double), dimension(:), intent(out), allocatable | values, | ||
integer(c_int), dimension(:), intent(out), allocatable | indices, | ||
type(error), intent(out), optional | err | ||
) |
Returns a copy of the specified global row in user-provided arrays.
[in] | globalrow | (In) Global row to extract. |
[out] | values | (Out) Extracted values for this row. |
[out] | indices | (Out) Extracted global column indices for the corresponding values. |
[out] | err | Integer error code, set to 0 if successful, non-zero if global row is not owned by calling process or if the number of entries in this row exceed the Length parameter. |
subroutine fepetra_crsmatrix::fillcomplete | ( | class(epetra_crsmatrix), intent(in) | this, |
logical, intent(in), optional | OptimizeDataStorage, | ||
type(error), intent(out), optional | err | ||
) |
Signal that data entry is complete. Perform transformations to local index space. This version of FillComplete assumes that the domain and range distributions are identical to the matrix row distributions.
[in] | optimizedatastorage | (In) If true, storage will be packed for optimal performance. |
[out] | err | error code, 0 if successful. Returns a positive warning code of 3 if the matrix is rectangular (meaning that the other overloading of FillComplete should have been called, with differen domain-map and range-map specified). |
subroutine fepetra_crsmatrix::fillcomplete | ( | class(epetra_crsmatrix), intent(in) | this, |
class(epetra_map), intent(in) | DomainMap, | ||
class(epetra_map), intent(in) | RangeMap, | ||
logical, intent(in), optional | OptimizeDataStorage, | ||
type(error), intent(out), optional | err | ||
) |
Signal that data entry is complete. Perform transformations to local index space. This version of FillComplete requires the explicit specification of the domain and range distribution maps. These maps are used for importing and exporting vector and multi-vector elements that are needed for distributed matrix computations.
[in] | domainmap | (In) Map that describes the distribution of vector and multi-vectors in the matrix domain. |
[in] | rangemap | (In) Map that describes the distribution of vector and multi-vectors in the matrix range. |
[in] | optimizedatastorage | (In) If true, storage will be packed for optimal performance.By default storage will be optimized. If you cannot tolerate the increased temporary memory use, should set this value to false. |
[out] | err | (Out) error code, 0 if successful. positive warning code of 2 if it is detected that the matrix-graph got out of sync since this matrix was constructed (for instance if graph. FillComplete() was called by another matrix that shares the graph) |
subroutine fepetra_crsmatrix::insertglobalvalues | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | GlobalRow, | ||
real(c_double), dimension(:), intent(in) | values, | ||
integer(c_int), dimension(:), intent(in) | indices, | ||
type(error), intent(out), optional | err | ||
) |
Insert a list of elements in a given global row of the matrix. This method is used to construct a matrix for the first time. It cannot be used if the matrix structure has already been fixed (via a call to FillComplete()).
[in] | globalrow | (In) Row number (in global coordinates) to put elements. |
[in] | values | Values to enter. |
[in] | indices | (In) Global column indices corresponding to values. |
[out] | err | (Out) Integer error code, set to 0 if successful. |
integer(c_int) function fepetra_crsmatrix::maxnumentries | ( | class(epetra_crsmatrix), intent(in) | this | ) |
Returns the maximum number of nonzero entries across all rows on this processor.
subroutine fepetra_crsmatrix::multiply_multivector | ( | class(epetra_crsmatrix), intent(in) | this, |
logical, intent(in) | TransA, | ||
class(epetra_multivector), intent(in) | x, | ||
class(epetra_multivector), intent(inout) | y, | ||
type(error), intent(inout), optional | err | ||
) |
Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_MultiVector X in Y.
[in] | transa | (In) If true, multiply by the transpose of matrix, otherwise just use matrix. |
[in] | x | (In) An Epetra_MultiVector of dimension NumVectors to multiply with matrix. |
[in,out] | y | (Out) An Epetra_MultiVector of dimension NumVectorscontaining result. |
[in,out] | err | (Out) Integer error code, set to 0 if successful. |
subroutine fepetra_crsmatrix::multiply_vector | ( | class(epetra_crsmatrix), intent(in) | this, |
logical, intent(in) | TransA, | ||
class(epetra_vector), intent(in) | x, | ||
class(epetra_vector), intent(inout) | y, | ||
type(error), intent(inout), optional | err | ||
) |
Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_Vector x in y.
[in] | transa | (In) If true, multiply by the transpose of matrix, otherwise just use matrix. |
[in] | x | (In) An Epetra_Vector to multiply by. |
[in,out] | y | (Out) An Epetra_Vector containing result. |
[in,out] | err | (Out) Integer error code, set to 0 if successful. |
logical function fepetra_crsmatrix::myglobalrow | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | GID | ||
) |
Returns true of GID is owned by the calling processor, otherwise it returns false.
integer(c_int) function fepetra_crsmatrix::numglobalentries | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | GlobalRow | ||
) |
Returns the current number of nonzero entries in specified global row on this processor.
integer(c_int) function fepetra_crsmatrix::nummyentries | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | MyRow | ||
) |
Returns the current number of nonzero entries in specified local row on this processor.
integer(c_int) function fepetra_crsmatrix::nummyrowentries | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | MyRow | ||
) |
Return the current number of values stored for the specified local row.
[in] | myrow | (In) Local row. |
subroutine fepetra_crsmatrix::putscalar | ( | class(epetra_crsmatrix), intent(in) | this, |
real(c_double), intent(in) | scalar, | ||
type(error), intent(out), optional | err | ||
) |
Initialize all values in the matrix with constant value.
[in] | scalar | (In) Value to use. |
[out] | err | (Out) Integer error code, set to 0 if successful. |
subroutine fepetra_crsmatrix::replaceglobalvalues | ( | class(epetra_crsmatrix), intent(in) | this, |
integer(c_int), intent(in) | GlobalRow, | ||
real(c_double), dimension(:) | values, | ||
integer(c_int), dimension(:) | indices, | ||
type(error), intent(out), optional | err | ||
) |
Replace specified existing values with this list of entries for a given global row of the matrix.
[in] | globalrow | (In) Row number (in global coordinates) to put elements. |
values | (In) Values to enter. | |
indices | (In) Global column indices corresponding to values. | |
[out] | err | Integer error code, set to 0 if successful. |
type(epetra_map) function fepetra_crsmatrix::rowmap | ( | class(epetra_crsmatrix), intent(in) | this | ) |
Returns the Epetra_Map object associated with the rows of this matrix.
type(epetra_map) function fepetra_crsmatrix::rowmatrixrowmap | ( | class(epetra_crsmatrix), intent(in) | this | ) |
Returns the Epetra_Map object associated with the rows of this matrix.
subroutine fepetra_crsmatrix::scale | ( | class(epetra_crsmatrix), intent(in) | this, |
real(c_double), intent(in) | ScalarConstant, | ||
type(error), intent(out), optional | err | ||
) |
Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A).
[in] | scalarconstant | (In) Value to use. |
[out] | err | (Out) Integer error code, set to 0 if successful. |