LifeV
MatrixEpetra< DataType > Class Template Reference

MatrixEpetra - The Epetra Matrix format Wrapper. More...

#include <MatrixEpetra.hpp>

+ Inheritance diagram for MatrixEpetra< DataType >:
+ Collaboration diagram for MatrixEpetra< DataType >:

Private Attributes

std::shared_ptr< MapEpetraM_map
 
std::shared_ptr< const MapEpetraM_domainMap
 
std::shared_ptr< const MapEpetraM_rangeMap
 Shared pointer on the range MapEpetra. More...
 
matrix_ptrtype M_epetraCrs
 

Public Types

typedef Epetra_FECrsMatrix matrix_type
 
typedef std::shared_ptr< matrix_typematrix_ptrtype
 
typedef VectorEpetra vector_type
 
typedef std::shared_ptr< vector_typevectorPtr_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

MatrixEpetraoperator+= (const MatrixEpetra &matrix)
 Addition operator. More...
 
MatrixEpetraoperator-= (const MatrixEpetra &matrix)
 Subtraction operator. More...
 
MatrixEpetraoperator= (const MatrixEpetra &matrix)
 Assignment operator. More...
 
vector_type operator* (const vector_type &vector) const
 Matrix-Vector multiplication. More...
 
MatrixEpetraoperator*= (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:

\[ C += A \otimes B + \]

, where

\[ C \]

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_ptrtypematrixPtr ()
 Return the shared_pointer of the Epetra_FECrsMatrix. More...
 
std::shared_ptr< MapEpetramapPtr ()
 Return the shared_pointer of the Epetra_Map. More...
 
const matrix_ptrtypematrixPtr () 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 MapEpetramap () const
 Return the row MapEpetra of the MatrixEpetra used in the assembling. More...
 
const MapEpetradomainMap () const
 Return the domain MapEpetra of the MatrixEpetra. More...
 
const boost::shared_ptr< const MapEpetra > & domainMapPtr () const
 
const MapEpetrarangeMap () 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...
 

Detailed Description

template<typename DataType>
class LifeV::MatrixEpetra< DataType >

MatrixEpetra - The Epetra Matrix format Wrapper.

Author
Gilles Fourestey gille.nosp@m.s.fo.nosp@m.urest.nosp@m.ey@e.nosp@m.pfl.c.nosp@m.h
Simone Deparis simon.nosp@m.e.de.nosp@m.paris.nosp@m.@epf.nosp@m.l.ch
Gwenol Grandperrin gweno.nosp@m.l.gr.nosp@m.andpe.nosp@m.rrin.nosp@m.@epfl.nosp@m..ch

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.

Member Typedef Documentation

◆ matrix_type

typedef Epetra_FECrsMatrix matrix_type

Definition at line 84 of file MatrixEpetra.hpp.

◆ matrix_ptrtype

typedef std::shared_ptr<matrix_type> matrix_ptrtype

Definition at line 85 of file MatrixEpetra.hpp.

◆ vector_type

Definition at line 86 of file MatrixEpetra.hpp.

◆ vectorPtr_Type

typedef std::shared_ptr<vector_type> vectorPtr_Type

Definition at line 87 of file MatrixEpetra.hpp.

Constructor & Destructor Documentation

◆ MatrixEpetra() [1/7]

MatrixEpetra ( const MapEpetra map,
const Epetra_CrsGraph &  graph,
bool  ignoreNonLocalValues = false 
)

Constructor from a graph.

Parameters
mapRow map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...)
graphA sparse compressed row graph.

Definition at line 655 of file MatrixEpetra.hpp.

◆ MatrixEpetra() [2/7]

MatrixEpetra ( const MapEpetra map,
const Epetra_FECrsGraph &  graph,
bool  ignoreNonLocalValues = false 
)

Constructor from a FE graph.

Parameters
mapRow map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...)
graphA sparse compressed row FE graph.

Definition at line 663 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ MatrixEpetra() [3/7]

MatrixEpetra ( const MapEpetra map,
Int  numEntries = 50,
bool  ignoreNonLocalValues = false 
)

Constructor for square and rectangular matrices.

Parameters
mapRow map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...)
numEntriesThe average number of entries for each row.

Definition at line 671 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ MatrixEpetra() [4/7]

MatrixEpetra ( const MapEpetra map,
Int numEntriesPerRow,
bool  ignoreNonLocalValues = false 
)

Constructor for square and rectangular matrices, knowing the number of entries per row.

Parameters
mapRow map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...)
numEntriesPerRowContains the number of entries for each row.

◆ MatrixEpetra() [5/7]

MatrixEpetra ( const MatrixEpetra< DataType > &  matrix)

Copy Constructor.

Parameters
matrixMatrix used to create the new occurence

Definition at line 687 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ MatrixEpetra() [6/7]

MatrixEpetra ( const MatrixEpetra< DataType > &  matrix,
const UInt  reduceToProc 
)

Copies the matrix to a matrix which resides only on the processor.

Parameters
matrixMatrix where the content of the matrix should be stored
reduceToProcProcessor where the matrix should be stored

Definition at line 697 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ MatrixEpetra() [7/7]

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

Parameters
mapRow map. The column map will be defined in MatrixEpetra<DataType>::GlobalAssemble(...,...)
crsMatrixPtrPointer on a Epetra_FECrsMatrix of Trilinos

Definition at line 723 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ ~MatrixEpetra()

virtual ~MatrixEpetra ( )
inlinevirtual

Destructor.

Definition at line 154 of file MatrixEpetra.hpp.

Member Function Documentation

◆ LIFEV_DEPRECATED()

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

Parameters
crsMatrixPtrPointer on a Epetra_FECrsMatrix of Trilinos

◆ operator+=()

MatrixEpetra< DataType > & operator+= ( const MatrixEpetra< DataType > &  matrix)

Addition operator.

Parameters
matrixmatrix to be added to the current matrix

Definition at line 738 of file MatrixEpetra.hpp.

◆ operator-=()

MatrixEpetra< DataType > & operator-= ( const MatrixEpetra< DataType > &  matrix)

Subtraction operator.

Parameters
matrixmatrix to be subtracted to the current matrix

Definition at line 747 of file MatrixEpetra.hpp.

◆ operator=()

MatrixEpetra< DataType > & operator= ( const MatrixEpetra< DataType > &  matrix)

Assignment operator.

Parameters
matrixmatrix to be assigned to the current matrix

Definition at line 755 of file MatrixEpetra.hpp.

◆ operator*() [1/2]

MatrixEpetra< DataType >::vector_type operator* ( const vector_type vector) const

Matrix-Vector multiplication.

Parameters
vectorVector to be multiplied by the current matrix

Definition at line 767 of file MatrixEpetra.hpp.

◆ operator*=()

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)

Parameters
scalarValue for the multiplication

Definition at line 785 of file MatrixEpetra.hpp.

◆ operator*() [2/2]

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)

Parameters
scalarValue for the multiplication

Definition at line 793 of file MatrixEpetra.hpp.

◆ openCrsMatrix()

void openCrsMatrix ( )

If the matrix has been filled, this function will reopen the Matrix.

Definition at line 803 of file MatrixEpetra.hpp.

◆ removeZeros()

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.

◆ swapCrsMatrix() [1/2]

void swapCrsMatrix ( matrix_ptrtype matrixPtr)

Swap the given shared pointer with the one of the matrix.

Parameters
matrixPtrpointer on the matrix

Definition at line 865 of file MatrixEpetra.hpp.

◆ swapCrsMatrix() [2/2]

void swapCrsMatrix ( MatrixEpetra< DataType > &  matrix)

Swap the matrix with the one given as argument.

Parameters
matrixmatrix which is swapped

Definition at line 872 of file MatrixEpetra.hpp.

◆ multiply() [1/2]

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.

Parameters
transposeCurrentIf true, it transposes the MatrixEpetra
matrixMatrix that multiply the MatrixEpetra
transposeMatrixif true, it transposes the matrix "matrix"
resultMatrix to store the result
callFillCompleteOnResultIf true, the matrix "result" will be filled (i.e. closed) after the multiplication

Definition at line 878 of file MatrixEpetra.hpp.

◆ multiply() [2/2]

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.

Parameters
transposeCurrentIf true, it transposes the MatrixEpetra
vector1Vector that will be multiply by the MatrixEpetra
vector2Vector that will store the result

Definition at line 913 of file MatrixEpetra.hpp.

◆ addDyadicProduct()

void addDyadicProduct ( const vector_type uniqueVector1,
const vector_type uniqueVector2 
)

Add to the matrix a dyadic product:

\[ C += A \otimes B + \]

, where

\[ C \]

is the matrix.

NOTE: This method has been tested only for square matrices and unique vectors

Parameters
vector1unique vector A
vector2unique vector B

Definition at line 927 of file MatrixEpetra.hpp.

◆ add()

void add ( const DataType  scalar,
const MatrixEpetra< DataType > &  matrix 
)

Add a multiple of a given matrix: *this += scalar*matrix.

Parameters
scalarScalar which multiplies the matrix
matrixMatrix to be added

Definition at line 962 of file MatrixEpetra.hpp.

◆ transpose()

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.

◆ diagonalize() [1/4]

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.

Parameters
rVecVector of the Id that should be set to "coefficient"
coefficientValue to be set on the diagonal
offsetOffset used for the indices

Definition at line 984 of file MatrixEpetra.hpp.

◆ diagonalize() [2/4]

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.

Parameters
entryIndexIndex of the row where we set the "coefficient"
coefficientValue to be set on the diagonal
offsetOffset used for the indices

Definition at line 1035 of file MatrixEpetra.hpp.

◆ diagonalize() [3/4]

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

Parameters
rVecvector of rows
coefficientValue to set entry (r,r) at
rhsRight hand side Vector of the system to be adapted accordingly
datumVectorvector of values to constrain entry r of the solution at
offsetOffset used for the indices

Definition at line 1075 of file MatrixEpetra.hpp.

◆ diagonalize() [4/4]

void diagonalize ( UInt const  row,
DataType const  coefficient,
vector_type rhs,
DataType  datum,
UInt  offset = 0 
)

apply constraint on row "row"

Parameters
rowrow number
coefficientvalue to set entry (row,row) at
rhsRight hand side vector of the system to be adapted accordingly
datumVectorvalue to constrain entry r of the solution at
offsetOffset used for the indices

Definition at line 1272 of file MatrixEpetra.hpp.

◆ matrixMarket()

void matrixMarket ( std::string const &  fileName,
const bool  headers = true 
)

Save the matrix into a MatrixMarket (.mtx) file.

Parameters
filenamefile where the matrix will be saved
headersboolean to write the MM headers or not

Definition at line 1330 of file MatrixEpetra.hpp.

◆ spy()

void spy ( std::string const &  fileName)

Save the matrix into a Matlab (.m) file.

Parameters
filenamefile where the matrix will be saved

Definition at line 1347 of file MatrixEpetra.hpp.

◆ exportToHDF5()

void exportToHDF5 ( std::string const &  fileName,
std::string const &  matrixName = "matrix",
bool const &  truncate = true 
)

Save the matrix into a HDF5 (.h5) file.

Parameters
fileNameName of the file where the matrix will be saved, without extension (.h5)
matrixNameName of the matrix in the HDF5 file
truncateTrue 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.

◆ importFromHDF5()

void importFromHDF5 ( std::string const &  fileName,
std::string const &  matrixName = "matrix" 
)

Read a matrix from a HDF5 (.h5) file.

Parameters
fileNameName of the file where the matrix will be saved, without extension (.h5)
matrixNameName of the matrix in the HDF5 file

Definition at line 1394 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ showMe()

void showMe ( std::ostream &  output = std::cout) const

Print the contents of the matrix.

Parameters
outputStream where the informations must be printed

Definition at line 1422 of file MatrixEpetra.hpp.

◆ globalAssemble() [1/2]

Int globalAssemble ( )

Global assemble of a square matrix with default domain and range map.

Definition at line 1428 of file MatrixEpetra.hpp.

◆ globalAssemble() [2/2]

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.

◆ fillComplete()

Int fillComplete ( )

Fill complete of a square matrix with default domain and range map.

Definition at line 1440 of file MatrixEpetra.hpp.

◆ insertValueDiagonal() [1/2]

void insertValueDiagonal ( const DataType  entry,
const MapEpetra Map,
const UInt  offset = 0 
)

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.

Parameters
entryThe entry that is inserted in the diagonal
MapThe MapEpetra
offsetAn offset for the insertion of the diagonal entries

Definition at line 1468 of file MatrixEpetra.hpp.

◆ insertValueDiagonal() [2/2]

void insertValueDiagonal ( const DataType &  value,
Int  from = -1,
Int  to = -2 
)

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;

  1. If from > to, process all diagonal entries entries
  2. If from = to, do nothing

This methods works only if matrix is not closed.

Parameters
valueValue to be inserted on the diagonal
fromStarting row
toEnding row

Definition at line 1477 of file MatrixEpetra.hpp.

◆ insertOneDiagonal()

void insertOneDiagonal ( Int  from = -1,
Int  to = -2 
)

insert ones into the diagonal to ensure the matrix' graph has a entry there

Definition at line 1519 of file MatrixEpetra.hpp.

◆ insertZeroDiagonal()

void insertZeroDiagonal ( Int  from = -1,
Int  to = -2 
)

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;

  1. If from > to, process all diagonal entries entries
  2. If from = to, do nothing

This methods works only if matrix is not closed.

Parameters
fromStarting row
toEnding row

Definition at line 1525 of file MatrixEpetra.hpp.

+ Here is the caller graph for this function:

◆ norm1()

Real norm1 ( ) const

Compute the norm 1 of the global matrix.

Returns
norm 1

Definition at line 1531 of file MatrixEpetra.hpp.

◆ normInf()

Real normInf ( ) const

Compute the norm inf of the global matrix.

Returns
norm inf

Definition at line 1537 of file MatrixEpetra.hpp.

◆ normFrobenius()

Real normFrobenius ( ) const

Compute the frobenius norm of the global matrix.

Returns
norm Frobenius

Definition at line 1543 of file MatrixEpetra.hpp.

◆ zero()

void zero ( )
inline

set zero in all the matrix entries

Definition at line 466 of file MatrixEpetra.hpp.

◆ setCoefficients()

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.

Parameters
numRowsNumber of rows into the list given in "localValues"
numColumnsNumber of columns into the list given in "localValues"
rowIndicesList of row indices
columnIndicesList of column indices
localValues2D array containing the coefficient related to "rowIndices" and "columnIndices"
formatFormat of the matrix (Epetra_FECrsMatrix::COLUMN_MAJOR or Epetra_FECrsMatrix::ROW_MAJOR)

Definition at line 1553 of file MatrixEpetra.hpp.

◆ setCoefficient()

void setCoefficient ( UInt  row,
UInt  column,
DataType  localValue 
)

Set a coefficient of the matrix.

Parameters
rowRow index of the coefficient
columncolumn index of the coefficient

Definition at line 1569 of file MatrixEpetra.hpp.

◆ addToCoefficients()

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.

Parameters
numRowsNumber of rows into the list given in "localValues"
numColumnsNumber of columns into the list given in "localValues"
rowIndicesList of row indices
columnIndicesList of column indices
localValues2D array containing the coefficient related to "rowIndices" and "columnIndices"
formatFormat of the matrix (Epetra_FECrsMatrix::COLUMN_MAJOR or Epetra_FECrsMatrix::ROW_MAJOR)

Definition at line 1604 of file MatrixEpetra.hpp.

◆ sumIntoCoefficients()

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.

Parameters
numRowsNumber of rows into the list given in "localValues"
numColumnsNumber of columns into the list given in "localValues"
rowIndicesList of row indices
columnIndicesList of column indices
localValues2D array containing the coefficient related to "rowIndices" and "columnIndices"
formatFormat of the matrix (Epetra_FECrsMatrix::COLUMN_MAJOR or Epetra_FECrsMatrix::ROW_MAJOR)

Definition at line 1625 of file MatrixEpetra.hpp.

◆ addToCoefficient()

void addToCoefficient ( UInt const  row,
UInt const  column,
DataType const  localValue 
)

Add a value at a coefficient of the matrix.

Parameters
rowRow index of the value to be added
columnColumn index of the value to be added
localValueValue to be added to the coefficient

Definition at line 1584 of file MatrixEpetra.hpp.

◆ setUseTranspose()

void setUseTranspose ( const bool &  transpose = false)
inline

If set true, transpose of this operator will be applied.

Parameters
transposeflag to identify whether to transpose or not

Definition at line 533 of file MatrixEpetra.hpp.

◆ filled()

bool filled ( ) const
inline

Return the fill-complete status of the Epetra_FECrsMatrix.

Definition at line 545 of file MatrixEpetra.hpp.

◆ matrixPtr() [1/2]

matrix_ptrtype& matrixPtr ( )
inline

Return the shared_pointer of the Epetra_FECrsMatrix.

Definition at line 551 of file MatrixEpetra.hpp.

◆ mapPtr()

std::shared_ptr<MapEpetra> mapPtr ( )
inline

Return the shared_pointer of the Epetra_Map.

Definition at line 557 of file MatrixEpetra.hpp.

◆ matrixPtr() [2/2]

const matrix_ptrtype& matrixPtr ( ) const
inline

Return the const shared_pointer of the Epetra_FECrsMatrix.

Definition at line 563 of file MatrixEpetra.hpp.

◆ meanNumEntries()

Int meanNumEntries ( ) const

Return the mean number of entries in the matrix rows.

Definition at line 1650 of file MatrixEpetra.hpp.

◆ processorId()

Int processorId ( )

Return the Id of the processor.

Definition at line 1667 of file MatrixEpetra.hpp.

◆ map()

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.

◆ domainMap()

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.

◆ domainMapPtr()

const boost::shared_ptr< const MapEpetra > & domainMapPtr ( ) const

Definition at line 1687 of file MatrixEpetra.hpp.

◆ rangeMap()

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.

◆ rangeMapPtr()

const boost::shared_ptr< const MapEpetra > & rangeMapPtr ( ) const

Definition at line 1700 of file MatrixEpetra.hpp.

◆ restrict()

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.

Parameters
mapMapEpetra that contains the indices
matrix_outMatrix restricted

Definition at line 1747 of file MatrixEpetra.hpp.

Friends And Related Function Documentation

◆ RAP

MatrixEpetra<DType>* RAP ( const MatrixEpetra< DType > &  R,
const MatrixEpetra< DType > &  A,
const MatrixEpetra< DType > &  P 
)
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.

◆ PtAP

MatrixEpetra<DType>* PtAP ( const MatrixEpetra< DType > &  A,
const MatrixEpetra< DType > &  P 
)
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.

Field Documentation

◆ M_map

std::shared_ptr< MapEpetra > M_map
private

Definition at line 626 of file MatrixEpetra.hpp.

◆ M_domainMap

std::shared_ptr< const MapEpetra > M_domainMap
private

Definition at line 635 of file MatrixEpetra.hpp.

◆ M_rangeMap

std::shared_ptr< const MapEpetra > M_rangeMap
private

Shared pointer on the range MapEpetra.

Definition at line 644 of file MatrixEpetra.hpp.

◆ M_epetraCrs

matrix_ptrtype M_epetraCrs
private

Definition at line 647 of file MatrixEpetra.hpp.


The documentation for this class was generated from the following file: