Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
Tpetra::MatrixMatrix Namespace Reference

Distributed sparse matrix-matrix multiply and add. More...

Functions

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Multiply (const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool transposeA, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, bool transposeB, CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &C, bool call_FillComplete_on_result=true, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Sparse matrix-matrix multiply. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Add (const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool transposeA, Scalar scalarA, CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarB)
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
add (const Scalar &alpha, const bool transposeA, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const bool transposeB, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Compute the sparse matrix sum C = scalarA * Op(A) + scalarB * Op(B), where Op(X) is either X or its transpose. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void add (const Scalar &alpha, const bool transposeA, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const bool transposeB, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &C, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Compute the sparse matrix sum C = scalarA * Op(A) + scalarB * Op(B), where Op(X) is either X or its transpose. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Add (const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool transposeA, Scalar scalarA, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, bool transposeB, Scalar scalarB, Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > C)
 Compute the sparse matrix sum C = scalarA * Op(A) + scalarB * Op(B), where Op(X) is either X or its transpose. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Jacobi (Scalar omega, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Dinv, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &C, bool call_FillComplete_on_result=true, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 

Detailed Description

Distributed sparse matrix-matrix multiply and add.

This namespace includes functions for computing the sum or product of two distributed sparse matrices, each of which is represented as a Tpetra::CrsMatrix.

Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MatrixMatrix::Multiply ( const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
bool  transposeA,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
bool  transposeB,
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  C,
bool  call_FillComplete_on_result = true,
const std::string &  label = std::string(),
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Sparse matrix-matrix multiply.

Given CrsMatrix instances A and B, compute the product C = A*B, overwriting an existing CrsMatrix instance C with the result.

Precondition
All three matrices A, B, and C must have uniquely owned row Maps.
On input, C must have the same row Map as A.
A and B must be fill complete.
If C has a range Map on input, then A and C must have the same range Maps.
If C has a domain Map on input, then B and C must have the same domain Maps.

For the latter two preconditions, recall that a matrix does not have a domain or range Map unless fillComplete has been called on it at least once.

Parameters
A[in] fill-complete sparse matrix.
transposeA[in] Whether to use transpose of matrix A.
B[in] fill-complete sparse matrix.
transposeB[in] Whether to use transpose of matrix B.
C[in/out] On entry to this method, if C is fill complete, then C's graph must have the correct structure, that is, its structure must equal the structure of A*B. On exit, C will be fill complete, unless the last argument to this function is false.
call_FillComplete_on_result[in] Optional argument; defaults to true. If false, C will not be fill complete on output.

Definition at line 100 of file TpetraExt_MatrixMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MatrixMatrix::Add ( const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
bool  transposeA,
Scalar  scalarA,
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
Scalar  scalarB 
)
Given CrsMatrix objects A and B, form the sum B = a*A + b*B

Currently not functional.

Parameters
AInput, must already have had 'FillComplete()' called.
transposeAInput, whether to use transpose of matrix A.
scalarAInput, scalar multiplier for matrix A.
BResult. On entry to this method, it doesn't matter whether FillComplete() has already been called on B or not. If it has, then B's graph must already contain all nonzero locations that will be produced when forming the sum.
scalarBInput, scalar multiplier for matrix B.

Definition at line 468 of file TpetraExt_MatrixMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Tpetra::MatrixMatrix::add ( const Scalar &  alpha,
const bool  transposeA,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Scalar &  beta,
const bool  transposeB,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Compute the sparse matrix sum C = scalarA * Op(A) + scalarB * Op(B), where Op(X) is either X or its transpose.

This version of sparse matrix-matrix add returns a new CrsMatrix instance, rather than using an existing instance for the result. The returned matrix is always fill complete, with the given domain and range Maps. It is correct (though less efficient) for A and B to have different row Maps; the returned matrix will have the same row Map as the row Map of B.

Precondition
A and B must both be fillComplete and have matching domain and range Maps.
Parameters
scalarA[in] Scalar multiplier for A in the sum.
transposeA[in] If true, use the transpose of A.
A[in] The first input matrix.
scalarB[in] Scalar multiplier for B in the sum.
transposeB[in] If true, use the transpose of B.
B[in] The second input matrix.
domainMap[in] Domain Map of C (on output). If null or not provided, this defaults to the domain map of Op(B).
rangeMap[in] Range Map of C (on output). If null or not provided, this defaults to the range map of Op(B).
params[in/out] Same as the parameters of RowMatrix::add.

See the documentation of RowMatrix::add for a more detailed discussion of the optional parameters.

Definition at line 542 of file TpetraExt_MatrixMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MatrixMatrix::add ( const Scalar &  alpha,
const bool  transposeA,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const Scalar &  beta,
const bool  transposeB,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  C,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Compute the sparse matrix sum C = scalarA * Op(A) + scalarB * Op(B), where Op(X) is either X or its transpose.

This version of sparse matrix-matrix add returns a new CrsMatrix instance, rather than using an existing instance for the result. The returned matrix is always locally indexed and is fill complete by default (with the given domain and range Maps, or using those of Op(B) if null). If params->get("Call fillComplete") = false, then the resulting matrix will not be fill complete but it will have an immutable sparsity pattern (only scalar values can be changed). This is because C is allocated with the exact amount of storage to hold the sum.

It is correct (though less efficient) for A and B to have different row Maps; the returned matrix will have the same row Map as the row Map of B.

Precondition
A and B must both be fillComplete and have matching domain and range Maps.
Parameters
scalarA[in] Scalar multiplier for A in the sum.
transposeA[in] If true, use the transpose of A.
A[in] The first input matrix.
scalarB[in] Scalar multiplier for B in the sum.
transposeB[in] If true, use the transpose of B.
B[in] The second input matrix.
C[out] The result matrix, which we expect to be 'new' (no entries inserted) on input.
domainMap[in] Domain Map of C (on output). If null or not provided, this defaults to the domain map of B (or range map, if transposeB).
rangeMap[in] Range Map of C (on output). If null or not provided, this defaults to the range map of B (or domain map, if transposeB).
params[in/out] Same as the parameters of RowMatrix::add.

See the documentation of RowMatrix::add for a more detailed discussion of the optional parameters.

Definition at line 610 of file TpetraExt_MatrixMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MatrixMatrix::Add ( const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
bool  transposeA,
Scalar  scalarA,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
bool  transposeB,
Scalar  scalarB,
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > >  C 
)

Compute the sparse matrix sum C = scalarA * Op(A) + scalarB * Op(B), where Op(X) is either X or its transpose.

Precondition
Both input matrices A and B must be fill complete. That is, their fillComplete() method must have been called at least once, without an intervening call to resumeFill().
Parameters
A[in] The first input matrix.
transposeA[in] If true, use the transpose of A.
scalarA[in] Scalar multiplier for A in the sum.
B[in] The second input matrix.
transposeB[in] If true, use the transpose of B.
scalarB[in] Scalar multiplier for B in the sum.
C[in/out] On entry, C may be either null or a valid matrix. If C is null on input, this function will allocate a new CrsMatrix to contain the sum. If C is not null and is fill complete, then this function assumes that the sparsity pattern of the sum is fixed and compatible with the sparsity pattern of A + B. If C is not null and is not fill complete, then this function returns without calling fillComplete on C.
Warning
The case where C == null on input does not actually work. In order for it to work, we would need to change the interface of this function (for example, to pass in C as a (pointer or nonconst reference) to a Teuchos::RCP). Please use add() (which see) if you want matrix-matrix add to return a new instance of CrsMatrix.

Definition at line 907 of file TpetraExt_MatrixMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::MatrixMatrix::Jacobi ( Scalar  omega,
const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Dinv,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  A,
const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  B,
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  C,
bool  call_FillComplete_on_result = true,
const std::string &  label = std::string(),
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Given CrsMatrix objects A, B and C, and Vector Dinv, form the product C = (I-omega * Dinv A)*B In a parallel setting, A and B need not have matching distributions, but C needs to have the same row-map as A.

Parameters
omegaInput, scalar multiplier for Dinverse A
DinvInput, Vector representing a diagonal matrix, must match A->getRowMap()
AInput, must already have had 'fillComplete()' called.
BInput, must already have had 'fillComplete()' called.
CResult. On entry to this method, it doesn't matter whether fillComplete() has already been called on C or not. If it has, then C's graph must already contain all nonzero locations that will be produced when forming the product A*B. On exit, C.FillComplete() will have been called, unless the last argument to this function is specified to be false.
call_fillComplete_on_resultOptional argument, defaults to true. Power users may specify this argument to be false if they DON'T want this function to call C.fillComplete. (It is often useful to allow this function to call C.fillComplete, in cases where one or both of the input matrices are rectangular and it is not trivial to know which maps to use for the domain- and range-maps.)

Definition at line 293 of file TpetraExt_MatrixMatrix_def.hpp.