LifeV
|
MatrixEpetra - The Epetra Matrix format Wrapper. More...
#include <MatrixEpetra.hpp>
Private Attributes | |
std::shared_ptr< MapEpetra > | M_map |
std::shared_ptr< const MapEpetra > | M_domainMap |
std::shared_ptr< const MapEpetra > | M_rangeMap |
Shared pointer on the range MapEpetra. More... | |
matrix_ptrtype | M_epetraCrs |
Public Types | |
typedef Epetra_FECrsMatrix | matrix_type |
typedef std::shared_ptr< matrix_type > | matrix_ptrtype |
typedef VectorEpetra | vector_type |
typedef std::shared_ptr< vector_type > | vectorPtr_Type |
Constructors & Destructor | |
MatrixEpetra (const MapEpetra &map, const Epetra_CrsGraph &graph, bool ignoreNonLocalValues=false) | |
Constructor from a graph. More... | |
MatrixEpetra (const MapEpetra &map, const Epetra_FECrsGraph &graph, bool ignoreNonLocalValues=false) | |
Constructor from a FE graph. More... | |
MatrixEpetra (const MapEpetra &map, Int numEntries=50, bool ignoreNonLocalValues=false) | |
Constructor for square and rectangular matrices. More... | |
MatrixEpetra (const MapEpetra &map, Int *numEntriesPerRow, bool ignoreNonLocalValues=false) | |
Constructor for square and rectangular matrices, knowing the number of entries per row. More... | |
MatrixEpetra (const MatrixEpetra &matrix) | |
Copy Constructor. More... | |
MatrixEpetra (const MatrixEpetra &matrix, const UInt reduceToProc) | |
Copies the matrix to a matrix which resides only on the processor. More... | |
LIFEV_DEPRECATED (MatrixEpetra(matrix_ptrtype crsMatrixPtr)) | |
Constructs an MatrixEpetra view of an Epetra_FECrsMatrix. More... | |
MatrixEpetra (const MapEpetra &map, matrix_ptrtype crsMatrixPtr) | |
Constructs an MatrixEpetra view of an Epetra_FECrsMatrix. More... | |
virtual | ~MatrixEpetra () |
Destructor. More... | |
Operators | |
MatrixEpetra & | operator+= (const MatrixEpetra &matrix) |
Addition operator. More... | |
MatrixEpetra & | operator-= (const MatrixEpetra &matrix) |
Subtraction operator. More... | |
MatrixEpetra & | operator= (const MatrixEpetra &matrix) |
Assignment operator. More... | |
vector_type | operator* (const vector_type &vector) const |
Matrix-Vector multiplication. More... | |
MatrixEpetra & | operator*= (const DataType scalar) |
Multiplication operator. More... | |
MatrixEpetra | operator* (const DataType scalar) const |
Multiplication operator. More... | |
Methods | |
void | openCrsMatrix () |
If the matrix has been filled, this function will reopen the Matrix. More... | |
void | removeZeros () |
This function removes all the zeros in the matrix and add zero on the diagonal. More... | |
void | swapCrsMatrix (matrix_ptrtype &matrixPtr) |
Swap the given shared pointer with the one of the matrix. More... | |
void | swapCrsMatrix (MatrixEpetra< DataType > &matrix) |
Swap the matrix with the one given as argument. More... | |
Int | multiply (bool transposeCurrent, const MatrixEpetra< DataType > &matrix, bool transposeMatrix, MatrixEpetra< DataType > &result, bool callFillCompleteOnResult=true) const |
Multiply the MatrixEpetra by the first given matrix and put the result in the second given matrix. More... | |
Int | multiply (bool transposeCurrent, const vector_type &vector1, vector_type &vector2) const |
Multiply the first VectorEpetra given as a parameter by the MatrixEpetra and put the result into the second given VectorEpetra. More... | |
void | addDyadicProduct (const vector_type &uniqueVector1, const vector_type &uniqueVector2) |
Add to the matrix a dyadic product:
, where
is the matrix. More... | |
void | add (const DataType scalar, const MatrixEpetra &matrix) |
Add a multiple of a given matrix: *this += scalar*matrix. More... | |
std::shared_ptr< MatrixEpetra< DataType > > | transpose () |
Returns a pointer to a new matrix which contains the transpose of the current matrix. More... | |
void | diagonalize (std::vector< UInt > rVec, DataType const coefficient, UInt offset=0) |
Set entries (rVec(i),rVec(i)) to coefficient and the rest of the row entries to zero. More... | |
void | diagonalize (UInt const entryIndex, DataType const coefficient, UInt offset=0) |
Set entry (entryIndex,entryIndex) to coefficient and the rest of the row to zero. More... | |
void | diagonalize (std::vector< UInt > rVec, DataType const coefficient, vector_type &rhs, std::vector< DataType > datumVector, UInt offset=0) |
apply constraint on all rows rVec More... | |
void | diagonalize (UInt const row, DataType const coefficient, vector_type &rhs, DataType datum, UInt offset=0) |
apply constraint on row "row" More... | |
void | matrixMarket (std::string const &fileName, const bool headers=true) |
Save the matrix into a MatrixMarket (.mtx) file. More... | |
void | spy (std::string const &fileName) |
Save the matrix into a Matlab (.m) file. More... | |
void | exportToHDF5 (std::string const &fileName, std::string const &matrixName="matrix", bool const &truncate=true) |
Save the matrix into a HDF5 (.h5) file. More... | |
void | importFromHDF5 (std::string const &fileName, std::string const &matrixName="matrix") |
Read a matrix from a HDF5 (.h5) file. More... | |
void | showMe (std::ostream &output=std::cout) const |
Print the contents of the matrix. More... | |
Int | globalAssemble () |
Global assemble of a square matrix with default domain and range map. More... | |
Int | globalAssemble (const std::shared_ptr< const MapEpetra > &domainMap, const std::shared_ptr< const MapEpetra > &rangeMap) |
Global assemble for rectangular matrices. More... | |
Int | fillComplete () |
Fill complete of a square matrix with default domain and range map. More... | |
void | insertValueDiagonal (const DataType entry, const MapEpetra &Map, const UInt offset=0) |
insert the given value into the diagonal More... | |
void | insertValueDiagonal (const DataType &value, Int from=-1, Int to=-2) |
insert the given value into the diagonal More... | |
void | insertOneDiagonal (Int from=-1, Int to=-2) |
insert ones into the diagonal to ensure the matrix' graph has a entry there More... | |
void | insertZeroDiagonal (Int from=-1, Int to=-2) |
insert zeros into the diagonal to ensure the matrix' graph has a entry there More... | |
Real | norm1 () const |
Compute the norm 1 of the global matrix. More... | |
Real | normInf () const |
Compute the norm inf of the global matrix. More... | |
Real | normFrobenius () const |
Compute the frobenius norm of the global matrix. More... | |
Set Methods | |
void | zero () |
set zero in all the matrix entries More... | |
void | setCoefficients (Int const numRows, Int const numColumns, std::vector< Int > const &rowIndices, std::vector< Int > const &columnIndices, DataType *const *const localValues, Int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
Set a set of values to the corresponding set of coefficient in the matrix. More... | |
void | setCoefficient (UInt row, UInt column, DataType localValue) |
Set a coefficient of the matrix. More... | |
void | addToCoefficients (Int const numRows, Int const numColumns, std::vector< Int > const &rowIndices, std::vector< Int > const &columnIndices, DataType *const *const localValues, Int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
Add a set of values to the corresponding set of coefficient in the matrix. More... | |
void | sumIntoCoefficients (Int const numRows, Int const numColumns, std::vector< Int > const &rowIndices, std::vector< Int > const &columnIndices, DataType *const *const localValues, Int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
Add a set of values to the corresponding set of coefficient in the closed matrix. More... | |
void | addToCoefficient (UInt const row, UInt const column, DataType const localValue) |
Add a value at a coefficient of the matrix. More... | |
void | setUseTranspose (const bool &transpose=false) |
If set true, transpose of this operator will be applied. More... | |
Get Methods | |
bool | filled () const |
Return the fill-complete status of the Epetra_FECrsMatrix. More... | |
matrix_ptrtype & | matrixPtr () |
Return the shared_pointer of the Epetra_FECrsMatrix. More... | |
std::shared_ptr< MapEpetra > | mapPtr () |
Return the shared_pointer of the Epetra_Map. More... | |
const matrix_ptrtype & | matrixPtr () const |
Return the const shared_pointer of the Epetra_FECrsMatrix. More... | |
Int | meanNumEntries () const |
Return the mean number of entries in the matrix rows. More... | |
Int | processorId () |
Return the Id of the processor. More... | |
const MapEpetra & | map () const |
Return the row MapEpetra of the MatrixEpetra used in the assembling. More... | |
const MapEpetra & | domainMap () const |
Return the domain MapEpetra of the MatrixEpetra. More... | |
const boost::shared_ptr< const MapEpetra > & | domainMapPtr () const |
const MapEpetra & | rangeMap () const |
Return the range MapEpetra of the MatrixEpetra. More... | |
const boost::shared_ptr< const MapEpetra > & | rangeMapPtr () const |
void | restrict (const std::shared_ptr< MapEpetra > &map, const std::shared_ptr< VectorEpetra > &numeration, const UInt &offset, std::shared_ptr< MatrixEpetra< DataType > > &matrix_out) |
Restrict the matrix to the dofs contained in the input map. More... | |
template<typename DType > | |
MatrixEpetra< DType > * | RAP (const MatrixEpetra< DType > &R, const MatrixEpetra< DType > &A, const MatrixEpetra< DType > &P) |
Friend Functions. More... | |
template<typename DType > | |
MatrixEpetra< DType > * | PtAP (const MatrixEpetra< DType > &A, const MatrixEpetra< DType > &P) |
PtAP matrix matrix multiplication result = Pt * A * P User is responsible to wrap the row pointer returned by this method with his favorite pointer. More... | |
MatrixEpetra - The Epetra Matrix format Wrapper.
The MatrixEpetra class provides a general interface for the Epetra_FECrsMatrix of Trilinos.
Visit http://trilinos.sandia.gov for more informations about Epetra_FECrsMatrix.
Definition at line 77 of file MatrixEpetra.hpp.
typedef Epetra_FECrsMatrix matrix_type |
Definition at line 84 of file MatrixEpetra.hpp.
typedef std::shared_ptr<matrix_type> matrix_ptrtype |
Definition at line 85 of file MatrixEpetra.hpp.
typedef VectorEpetra vector_type |
Definition at line 86 of file MatrixEpetra.hpp.
typedef std::shared_ptr<vector_type> vectorPtr_Type |
Definition at line 87 of file MatrixEpetra.hpp.
MatrixEpetra | ( | const MapEpetra & | map, |
const Epetra_CrsGraph & | graph, | ||
bool | ignoreNonLocalValues = false |
||
) |
Constructor from a graph.
map | Row map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...) |
graph | A sparse compressed row graph. |
Definition at line 655 of file MatrixEpetra.hpp.
MatrixEpetra | ( | const MapEpetra & | map, |
const Epetra_FECrsGraph & | graph, | ||
bool | ignoreNonLocalValues = false |
||
) |
Constructor from a FE graph.
map | Row map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...) |
graph | A sparse compressed row FE graph. |
Definition at line 663 of file MatrixEpetra.hpp.
MatrixEpetra | ( | const MapEpetra & | map, |
Int | numEntries = 50 , |
||
bool | ignoreNonLocalValues = false |
||
) |
Constructor for square and rectangular matrices.
map | Row map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...) |
numEntries | The average number of entries for each row. |
Definition at line 671 of file MatrixEpetra.hpp.
MatrixEpetra | ( | const MapEpetra & | map, |
Int * | numEntriesPerRow, | ||
bool | ignoreNonLocalValues = false |
||
) |
Constructor for square and rectangular matrices, knowing the number of entries per row.
map | Row map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...) |
numEntriesPerRow | Contains the number of entries for each row. |
MatrixEpetra | ( | const MatrixEpetra< DataType > & | matrix | ) |
Copy Constructor.
matrix | Matrix used to create the new occurence |
Definition at line 687 of file MatrixEpetra.hpp.
MatrixEpetra | ( | const MatrixEpetra< DataType > & | matrix, |
const UInt | reduceToProc | ||
) |
Copies the matrix to a matrix which resides only on the processor.
matrix | Matrix where the content of the matrix should be stored |
reduceToProc | Processor where the matrix should be stored |
Definition at line 697 of file MatrixEpetra.hpp.
MatrixEpetra | ( | const MapEpetra & | map, |
matrix_ptrtype | crsMatrixPtr | ||
) |
Constructs an MatrixEpetra view of an Epetra_FECrsMatrix.
This constructor can be used when we need to modify an Epetra_FECrsMatrix using a method of the class MatrixEpetra
map | Row map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...) |
crsMatrixPtr | Pointer on a Epetra_FECrsMatrix of Trilinos |
Definition at line 723 of file MatrixEpetra.hpp.
|
inlinevirtual |
Destructor.
Definition at line 154 of file MatrixEpetra.hpp.
LIFEV_DEPRECATED | ( | MatrixEpetra< DataType >(matrix_ptrtype crsMatrixPtr) | ) |
Constructs an MatrixEpetra view of an Epetra_FECrsMatrix.
This constructor can be used when we need to modify an Epetra_FECrsMatrix using a method of the class MatrixEpetra
crsMatrixPtr | Pointer on a Epetra_FECrsMatrix of Trilinos |
MatrixEpetra< DataType > & operator+= | ( | const MatrixEpetra< DataType > & | matrix | ) |
Addition operator.
matrix | matrix to be added to the current matrix |
Definition at line 738 of file MatrixEpetra.hpp.
MatrixEpetra< DataType > & operator-= | ( | const MatrixEpetra< DataType > & | matrix | ) |
Subtraction operator.
matrix | matrix to be subtracted to the current matrix |
Definition at line 747 of file MatrixEpetra.hpp.
MatrixEpetra< DataType > & operator= | ( | const MatrixEpetra< DataType > & | matrix | ) |
Assignment operator.
matrix | matrix to be assigned to the current matrix |
Definition at line 755 of file MatrixEpetra.hpp.
MatrixEpetra< DataType >::vector_type operator* | ( | const vector_type & | vector | ) | const |
Matrix-Vector multiplication.
vector | Vector to be multiplied by the current matrix |
Definition at line 767 of file MatrixEpetra.hpp.
MatrixEpetra< DataType > & operator*= | ( | const DataType | scalar | ) |
Multiplication operator.
Multiply by a scalar value the components of the current matrix. (modify the matrix of the class)
scalar | Value for the multiplication |
Definition at line 785 of file MatrixEpetra.hpp.
MatrixEpetra< DataType > operator* | ( | const DataType | scalar | ) | const |
Multiplication operator.
Multiply by a scalar value the components of the current matrix. (do not modify the matrix of the class)
scalar | Value for the multiplication |
Definition at line 793 of file MatrixEpetra.hpp.
void openCrsMatrix | ( | ) |
If the matrix has been filled, this function will reopen the Matrix.
Definition at line 803 of file MatrixEpetra.hpp.
void removeZeros | ( | ) |
This function removes all the zeros in the matrix and add zero on the diagonal.
The zeros added on the diagonal are used to impose boundary conditions
Definition at line 820 of file MatrixEpetra.hpp.
void swapCrsMatrix | ( | matrix_ptrtype & | matrixPtr | ) |
Swap the given shared pointer with the one of the matrix.
matrixPtr | pointer on the matrix |
Definition at line 865 of file MatrixEpetra.hpp.
void swapCrsMatrix | ( | MatrixEpetra< DataType > & | matrix | ) |
Swap the matrix with the one given as argument.
matrix | matrix which is swapped |
Definition at line 872 of file MatrixEpetra.hpp.
Int multiply | ( | bool | transposeCurrent, |
const MatrixEpetra< DataType > & | matrix, | ||
bool | transposeMatrix, | ||
MatrixEpetra< DataType > & | result, | ||
bool | callFillCompleteOnResult = true |
||
) | const |
Multiply the MatrixEpetra by the first given matrix and put the result in the second given matrix.
transposeCurrent | If true, it transposes the MatrixEpetra |
matrix | Matrix that multiply the MatrixEpetra |
transposeMatrix | if true, it transposes the matrix "matrix" |
result | Matrix to store the result |
callFillCompleteOnResult | If true, the matrix "result" will be filled (i.e. closed) after the multiplication |
Definition at line 878 of file MatrixEpetra.hpp.
Int multiply | ( | bool | transposeCurrent, |
const vector_type & | vector1, | ||
vector_type & | vector2 | ||
) | const |
Multiply the first VectorEpetra given as a parameter by the MatrixEpetra and put the result into the second given VectorEpetra.
transposeCurrent | If true, it transposes the MatrixEpetra |
vector1 | Vector that will be multiply by the MatrixEpetra |
vector2 | Vector that will store the result |
Definition at line 913 of file MatrixEpetra.hpp.
void addDyadicProduct | ( | const vector_type & | uniqueVector1, |
const vector_type & | uniqueVector2 | ||
) |
Add to the matrix a dyadic product:
, where
is the matrix.
NOTE: This method has been tested only for square matrices and unique vectors
vector1 | unique vector A |
vector2 | unique vector B |
Definition at line 927 of file MatrixEpetra.hpp.
void add | ( | const DataType | scalar, |
const MatrixEpetra< DataType > & | matrix | ||
) |
Add a multiple of a given matrix: *this += scalar*matrix.
scalar | Scalar which multiplies the matrix |
matrix | Matrix to be added |
Definition at line 962 of file MatrixEpetra.hpp.
std::shared_ptr< MatrixEpetra< DataType > > transpose | ( | ) |
Returns a pointer to a new matrix which contains the transpose of the current matrix.
Definition at line 968 of file MatrixEpetra.hpp.
Set entries (rVec(i),rVec(i)) to coefficient and the rest of the row entries to zero.
rVec | Vector of the Id that should be set to "coefficient" |
coefficient | Value to be set on the diagonal |
offset | Offset used for the indices |
Definition at line 984 of file MatrixEpetra.hpp.
Set entry (entryIndex,entryIndex) to coefficient and the rest of the row to zero.
entryIndex | Index of the row where we set the "coefficient" |
coefficient | Value to be set on the diagonal |
offset | Offset used for the indices |
Definition at line 1035 of file MatrixEpetra.hpp.
void diagonalize | ( | std::vector< UInt > | rVec, |
DataType const | coefficient, | ||
vector_type & | rhs, | ||
std::vector< DataType > | datumVector, | ||
UInt | offset = 0 |
||
) |
apply constraint on all rows rVec
rVec | vector of rows |
coefficient | Value to set entry (r,r) at |
rhs | Right hand side Vector of the system to be adapted accordingly |
datumVector | vector of values to constrain entry r of the solution at |
offset | Offset used for the indices |
Definition at line 1075 of file MatrixEpetra.hpp.
void diagonalize | ( | UInt const | row, |
DataType const | coefficient, | ||
vector_type & | rhs, | ||
DataType | datum, | ||
UInt | offset = 0 |
||
) |
apply constraint on row "row"
row | row number |
coefficient | value to set entry (row,row) at |
rhs | Right hand side vector of the system to be adapted accordingly |
datumVector | value to constrain entry r of the solution at |
offset | Offset used for the indices |
Definition at line 1272 of file MatrixEpetra.hpp.
void matrixMarket | ( | std::string const & | fileName, |
const bool | headers = true |
||
) |
Save the matrix into a MatrixMarket (.mtx) file.
filename | file where the matrix will be saved |
headers | boolean to write the MM headers or not |
Definition at line 1330 of file MatrixEpetra.hpp.
void spy | ( | std::string const & | fileName | ) |
Save the matrix into a Matlab (.m) file.
filename | file where the matrix will be saved |
Definition at line 1347 of file MatrixEpetra.hpp.
void exportToHDF5 | ( | std::string const & | fileName, |
std::string const & | matrixName = "matrix" , |
||
bool const & | truncate = true |
||
) |
Save the matrix into a HDF5 (.h5) file.
fileName | Name of the file where the matrix will be saved, without extension (.h5) |
matrixName | Name of the matrix in the HDF5 file |
truncate | True if the file has to be truncated; False if the file already exist and should not be truncated |
Definition at line 1363 of file MatrixEpetra.hpp.
void importFromHDF5 | ( | std::string const & | fileName, |
std::string const & | matrixName = "matrix" |
||
) |
Read a matrix from a HDF5 (.h5) file.
fileName | Name of the file where the matrix will be saved, without extension (.h5) |
matrixName | Name of the matrix in the HDF5 file |
Definition at line 1394 of file MatrixEpetra.hpp.
void showMe | ( | std::ostream & | output = std::cout | ) | const |
Print the contents of the matrix.
output | Stream where the informations must be printed |
Definition at line 1422 of file MatrixEpetra.hpp.
Int globalAssemble | ( | ) |
Global assemble of a square matrix with default domain and range map.
Definition at line 1428 of file MatrixEpetra.hpp.
Int globalAssemble | ( | const std::shared_ptr< const MapEpetra > & | domainMap, |
const std::shared_ptr< const MapEpetra > & | rangeMap | ||
) |
Global assemble for rectangular matrices.
Definition at line 1452 of file MatrixEpetra.hpp.
Int fillComplete | ( | ) |
Fill complete of a square matrix with default domain and range map.
Definition at line 1440 of file MatrixEpetra.hpp.
insert the given value into the diagonal
Pay intention that this will add values to the diagonal, so for later added values with set_mat_inc, the one will be added. Inserts Value on the local diagonal for diagonal elements specified by the input MapEpetra; This methods works only if matrix is not closed.
entry | The entry that is inserted in the diagonal |
Map | The MapEpetra |
offset | An offset for the insertion of the diagonal entries |
Definition at line 1468 of file MatrixEpetra.hpp.
insert the given value into the diagonal
Pay intention that this will add values to the diagonal, so for later added values with set_mat_inc, the one will be added
Inserts Value on the local diagonal for diagonal elements >= from and < to;
This methods works only if matrix is not closed.
value | Value to be inserted on the diagonal |
from | Starting row |
to | Ending row |
Definition at line 1477 of file MatrixEpetra.hpp.
insert ones into the diagonal to ensure the matrix' graph has a entry there
Definition at line 1519 of file MatrixEpetra.hpp.
insert zeros into the diagonal to ensure the matrix' graph has a entry there
This method does not remove non zero entries in the diagonal.
Inserts Value on the local diagonal for diagonal elements >= from and < to;
This methods works only if matrix is not closed.
from | Starting row |
to | Ending row |
Definition at line 1525 of file MatrixEpetra.hpp.
Real norm1 | ( | ) | const |
Compute the norm 1 of the global matrix.
Definition at line 1531 of file MatrixEpetra.hpp.
Real normInf | ( | ) | const |
Compute the norm inf of the global matrix.
Definition at line 1537 of file MatrixEpetra.hpp.
Real normFrobenius | ( | ) | const |
Compute the frobenius norm of the global matrix.
Definition at line 1543 of file MatrixEpetra.hpp.
|
inline |
set zero in all the matrix entries
Definition at line 466 of file MatrixEpetra.hpp.
void setCoefficients | ( | Int const | numRows, |
Int const | numColumns, | ||
std::vector< Int > const & | rowIndices, | ||
std::vector< Int > const & | columnIndices, | ||
DataType *const *const | localValues, | ||
Int | format = Epetra_FECrsMatrix::COLUMN_MAJOR |
||
) |
Set a set of values to the corresponding set of coefficient in the matrix.
numRows | Number of rows into the list given in "localValues" |
numColumns | Number of columns into the list given in "localValues" |
rowIndices | List of row indices |
columnIndices | List of column indices |
localValues | 2D array containing the coefficient related to "rowIndices" and "columnIndices" |
format | Format of the matrix (Epetra_FECrsMatrix::COLUMN_MAJOR or Epetra_FECrsMatrix::ROW_MAJOR) |
Definition at line 1553 of file MatrixEpetra.hpp.
Set a coefficient of the matrix.
row | Row index of the coefficient |
column | column index of the coefficient |
Definition at line 1569 of file MatrixEpetra.hpp.
void addToCoefficients | ( | Int const | numRows, |
Int const | numColumns, | ||
std::vector< Int > const & | rowIndices, | ||
std::vector< Int > const & | columnIndices, | ||
DataType *const *const | localValues, | ||
Int | format = Epetra_FECrsMatrix::COLUMN_MAJOR |
||
) |
Add a set of values to the corresponding set of coefficient in the matrix.
numRows | Number of rows into the list given in "localValues" |
numColumns | Number of columns into the list given in "localValues" |
rowIndices | List of row indices |
columnIndices | List of column indices |
localValues | 2D array containing the coefficient related to "rowIndices" and "columnIndices" |
format | Format of the matrix (Epetra_FECrsMatrix::COLUMN_MAJOR or Epetra_FECrsMatrix::ROW_MAJOR) |
Definition at line 1604 of file MatrixEpetra.hpp.
void sumIntoCoefficients | ( | Int const | numRows, |
Int const | numColumns, | ||
std::vector< Int > const & | rowIndices, | ||
std::vector< Int > const & | columnIndices, | ||
DataType *const *const | localValues, | ||
Int | format = Epetra_FECrsMatrix::COLUMN_MAJOR |
||
) |
Add a set of values to the corresponding set of coefficient in the closed matrix.
numRows | Number of rows into the list given in "localValues" |
numColumns | Number of columns into the list given in "localValues" |
rowIndices | List of row indices |
columnIndices | List of column indices |
localValues | 2D array containing the coefficient related to "rowIndices" and "columnIndices" |
format | Format of the matrix (Epetra_FECrsMatrix::COLUMN_MAJOR or Epetra_FECrsMatrix::ROW_MAJOR) |
Definition at line 1625 of file MatrixEpetra.hpp.
Add a value at a coefficient of the matrix.
row | Row index of the value to be added |
column | Column index of the value to be added |
localValue | Value to be added to the coefficient |
Definition at line 1584 of file MatrixEpetra.hpp.
|
inline |
If set true, transpose of this operator will be applied.
transpose | flag to identify whether to transpose or not |
Definition at line 533 of file MatrixEpetra.hpp.
|
inline |
Return the fill-complete status of the Epetra_FECrsMatrix.
Definition at line 545 of file MatrixEpetra.hpp.
|
inline |
Return the shared_pointer of the Epetra_FECrsMatrix.
Definition at line 551 of file MatrixEpetra.hpp.
|
inline |
Return the shared_pointer of the Epetra_Map.
Definition at line 557 of file MatrixEpetra.hpp.
|
inline |
Return the const shared_pointer of the Epetra_FECrsMatrix.
Definition at line 563 of file MatrixEpetra.hpp.
Int meanNumEntries | ( | ) | const |
Return the mean number of entries in the matrix rows.
Definition at line 1650 of file MatrixEpetra.hpp.
Int processorId | ( | ) |
Return the Id of the processor.
Definition at line 1667 of file MatrixEpetra.hpp.
const MapEpetra & map | ( | ) | const |
Return the row MapEpetra of the MatrixEpetra used in the assembling.
Note: This method should be call when MapEpetra is still open.
Definition at line 1673 of file MatrixEpetra.hpp.
const MapEpetra & domainMap | ( | ) | const |
Return the domain MapEpetra of the MatrixEpetra.
This function should be called only after MatrixEpetra<DataType>::GlobalAssemble(...) has been called. If this is an open matrix that M_domainMap is an invalid pointer
Definition at line 1680 of file MatrixEpetra.hpp.
const boost::shared_ptr< const MapEpetra > & domainMapPtr | ( | ) | const |
Definition at line 1687 of file MatrixEpetra.hpp.
const MapEpetra & rangeMap | ( | ) | const |
Return the range MapEpetra of the MatrixEpetra.
This function should be called only after MatrixEpetra<DataType>::GlobalAssemble(...) has been called. If this is an open matrix that M_domainMap is an invalid pointer
Definition at line 1693 of file MatrixEpetra.hpp.
const boost::shared_ptr< const MapEpetra > & rangeMapPtr | ( | ) | const |
Definition at line 1700 of file MatrixEpetra.hpp.
void restrict | ( | const std::shared_ptr< MapEpetra > & | map, |
const std::shared_ptr< VectorEpetra > & | numeration, | ||
const UInt & | offset, | ||
std::shared_ptr< MatrixEpetra< DataType > > & | matrix_out | ||
) |
Restrict the matrix to the dofs contained in the input map.
map | MapEpetra that contains the indices |
matrix_out | Matrix restricted |
Definition at line 1747 of file MatrixEpetra.hpp.
|
friend |
Friend Functions.
RAP matrix matrix multiplication result = R * A * P User is responsible to wrap the row pointer returned by this method with his favorite pointer
Definition at line 1706 of file MatrixEpetra.hpp.
|
friend |
PtAP matrix matrix multiplication result = Pt * A * P User is responsible to wrap the row pointer returned by this method with his favorite pointer.
Definition at line 1734 of file MatrixEpetra.hpp.
|
private |
Definition at line 626 of file MatrixEpetra.hpp.
|
private |
Definition at line 635 of file MatrixEpetra.hpp.
|
private |
Shared pointer on the range MapEpetra.
Definition at line 644 of file MatrixEpetra.hpp.
|
private |
Definition at line 647 of file MatrixEpetra.hpp.