Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps > Class Template Reference

VbrMatrix: Variable block row matrix. More...

#include <Tpetra_VbrMatrix_decl.hpp>

Inheritance diagram for Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >:
Inheritance graph
[legend]

List of all members.

Public Member Functions

Constructor/Destructor Methods

 VbrMatrix (const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &blkRowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile)
 Constructor specifying the row-map and the max number of (block) non-zeros for all rows.
virtual ~VbrMatrix ()
 Not Yet Implemented! Constructor specifying a pre-filled graph and block-maps for range and domain.
Advanced Matrix-vector multiplication method

template<class DomainScalar , class RangeScalar >
void multiply (const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp trans, RangeScalar alpha, RangeScalar beta) const
 Multiply this matrix by a MultiVector.
Operator Methods

const Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getDomainMap () const
 Returns the Map associated with the domain of this operator, which must be compatible with X.getMap().
const Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getRangeMap () const
 Returns the Map associated with the range of this operator, which must be compatible with Y.getMap().
void apply (const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp trans=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
 Computes the operator-multivector application.
bool hasTransposeApply () const
 Indicates whether this operator supports applying the adjoint operator.
Attribute Query Methods

const Teuchos::RCP< const
BlockMap< LocalOrdinal,
GlobalOrdinal, Node > > & 
getBlockRowMap () const
 Returns the block-row map.
const Teuchos::RCP< const
BlockMap< LocalOrdinal,
GlobalOrdinal, Node > > & 
getBlockColMap () const
 Returns the block-column map.
const Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getPointRowMap () const
 Returns the point-row map.
const Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getPointColMap () const
 Returns the point-column map.
bool isFillComplete () const
 Return true if fillComplete has been called, false otherwise.
Insertion Methods

void putScalar (Scalar s)
 Set the specified scalar throughout the matrix.
void setGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > &blockEntry)
 Copy the contents of the input block-entry into the matrix.
void sumIntoGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > &blockEntry)
 Add the contents of the input block-entry into the matrix.
void setGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView< const Scalar > &blockEntry)
 Copy the contents of the input block-entry into the matrix.
void sumIntoGlobalBlockEntry (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal blkRowSize, LocalOrdinal blkColSize, LocalOrdinal LDA, const Teuchos::ArrayView< const Scalar > &blockEntry)
 Add the contents of the input block-entry into the matrix.
Transformational Methods

void fillComplete (const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &blockDomainMap, const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &blockRangeMap)
void fillComplete ()
Extraction Methods

void getGlobalBlockEntryView (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< const Scalar > &blockEntry) const
 Returns a const read-only view of a block-entry.
void getGlobalBlockEntryViewNonConst (GlobalOrdinal globalBlockRow, GlobalOrdinal globalBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< Scalar > &blockEntry)
 Returns a non-const read-write view of a block-entry.
void getLocalBlockEntryView (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< const Scalar > &blockEntry) const
 Returns a const read-only view of a block-entry.
void getLocalBlockEntryViewNonConst (LocalOrdinal localBlockRow, LocalOrdinal localBlockCol, LocalOrdinal &numPtRows, LocalOrdinal &numPtCols, Teuchos::ArrayRCP< Scalar > &blockEntry)
 Returns a non-const read-write view of a block-entry.
Overridden from Teuchos::Describable

std::string description () const
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print the object with some verbosity level to a FancyOStream object.

Detailed Description

template<class Scalar, class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType, class LocalMatOps = typename Kokkos::DefaultKernels<Scalar,LocalOrdinal,Node>::BlockSparseOps>
class Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >

VbrMatrix: Variable block row matrix.

The VbrMatrix class has two significant 'states', distinguished by whether or not storage has been optimized (packed) or not.

When the matrix is in the non-optimized-storage state, internal data storage is in an un-packed, non-contiguous data-structure that allows for convenient insertion of data.

When the matrix is in the optimized-storage state, internal data is stored in contiguous (packed) arrays. When in this state, existing entries may be updated and replaced, but no new entries (indices and/or coefficients) may be inserted. In other words, the sparsity pattern or structure of the matrix may not be changed.

Use of the matrix as an Operator (performing matrix-vector multiplication) is only allowed when it is in the optimized-storage state.

VbrMatrix has two constructors, one which leaves the matrix in the optimized- storage stage, and another which leaves the matrix in the non-optimized-storage stage.

When the VbrMatrix is constructed in the non-optimized-storage state, and then filled using methods such as setGlobalBlockEntry etc., it can then be transformed to the optimized-storage state by calling the method fillComplete().

Once in the optimized-storage state, the VbrMatrix can not be returned to the non-optimized-storage state.

Definition at line 85 of file Tpetra_VbrMatrix_decl.hpp.


Constructor & Destructor Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::VbrMatrix ( const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > &  blkRowMap,
size_t  maxNumEntriesPerRow,
ProfileType  pftype = DynamicProfile 
) [inline]

Constructor specifying the row-map and the max number of (block) non-zeros for all rows.

After this constructor completes, the VbrMatrix is in the non-packed, non-optimized-storage, isFillComplete()==false state. Block-entries (rectangular, dense submatrices) may be inserted using class methods such as setGlobalBlockEntry(...), declared below.

Definition at line 138 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::~VbrMatrix (  )  [inline, virtual]

Not Yet Implemented! Constructor specifying a pre-filled graph and block-maps for range and domain.

Constructing a VbrMatrix with a pre-filled graph means that the matrix will start out in the optimized-storage, isFillComplete()==true state. The graph provided to this constructor must be already filled (If blkGraph->isFillComplete() != true, an exception is thrown.)

Entries in the input CrsGraph will correspond to block-entries in the VbrMatrix. In other words, the VbrMatrix will have a block-row corresponding to each row in the graph, and a block-entry corresponding to each column- index in the graph.

The block-maps provided for range and domain must be sized such that: blkDomainMap->getGlobalNumBlocks() == blkGraph->getDomainMap()->getGlobalNumElements(), blkRangeMap->getGlobalNumBlocks() == blkGraph->getRangeMap()->getGlobalNumElements(), blkDomainMap->getNodeNumBlocks() == blkGraph->getDomainMap()->getNodeNumElements(), blkRangeMap->getNodeNumBlocks() == blkGraph->getRangeMap()->getNodeNumElements(). If any of these conditions is not met, an exception is thrown. Destructor

Definition at line 175 of file Tpetra_VbrMatrix_def.hpp.


Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
template<class DomainScalar , class RangeScalar >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::multiply ( const MultiVector< DomainScalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
MultiVector< RangeScalar, LocalOrdinal, GlobalOrdinal, Node > &  Y,
Teuchos::ETransp  trans,
RangeScalar  alpha,
RangeScalar  beta 
) const [inline]

Multiply this matrix by a MultiVector.

X is required to be post-imported, i.e., described by the column map of the matrix. Y is required to be pre-exported, i.e., described by the row map of the matrix. See also the Operator::apply method which is implemented below.

Definition at line 207 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getDomainMap (  )  const [inline, virtual]

Returns the Map associated with the domain of this operator, which must be compatible with X.getMap().

Note that this is a point-entry map, not a block-map.

Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 182 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getRangeMap (  )  const [inline, virtual]

Returns the Map associated with the range of this operator, which must be compatible with Y.getMap().

Note that this is a point-entry map, not a block-map.

Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 190 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::apply ( const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  X,
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  Y,
Teuchos::ETransp  trans = Teuchos::NO_TRANS,
Scalar  alpha = Teuchos::ScalarTraits<Scalar>::one(),
Scalar  beta = Teuchos::ScalarTraits<Scalar>::zero() 
) const [inline, virtual]

Computes the operator-multivector application.

Loosely, performs $Y = \alpha \cdot A^{\textrm{trans}} \cdot X + \beta \cdot Y$. However, the details of operation vary according to the values of alpha and beta. Specifically

  • if beta == 0, apply() must overwrite Y, so that any values in Y (including NaNs) are ignored.
  • if alpha == 0, apply() may short-circuit the operator, so that any values in X (including NaNs) are ignored.

Implements Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 249 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
bool Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::hasTransposeApply (  )  const [inline, virtual]

Indicates whether this operator supports applying the adjoint operator.

Reimplemented from Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 274 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getBlockRowMap (  )  const [inline]

Returns the block-row map.

Definition at line 282 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
const Teuchos::RCP< const BlockMap< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getBlockColMap (  )  const [inline]

Returns the block-column map.

Definition at line 298 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getPointRowMap (  )  const [inline]

Returns the point-row map.

Definition at line 290 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getPointColMap (  )  const [inline]

Returns the point-column map.

Definition at line 306 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
bool Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::isFillComplete (  )  const [inline]

Return true if fillComplete has been called, false otherwise.

Definition at line 198 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::putScalar ( Scalar  s  )  [inline]

Set the specified scalar throughout the matrix.

This method may be called any time (before or after fillComplete()).

Definition at line 522 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setGlobalBlockEntry ( GlobalOrdinal  globalBlockRow,
GlobalOrdinal  globalBlockCol,
const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > &  blockEntry 
) [inline]

Copy the contents of the input block-entry into the matrix.

This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.

If the specified block-entry already exists in the matrix, it will be over-written (replaced) by the input block-entry.

This method may be called any time (before or after fillComplete()).

Definition at line 569 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoGlobalBlockEntry ( GlobalOrdinal  globalBlockRow,
GlobalOrdinal  globalBlockCol,
const Teuchos::SerialDenseMatrix< GlobalOrdinal, Scalar > &  blockEntry 
) [inline]

Add the contents of the input block-entry into the matrix.

This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.

If the specified block-entry already exists in the matrix, the contents of the input block-entry will be added to the values that are already present.

This method may be called any time (before or after fillComplete()).

Definition at line 585 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::setGlobalBlockEntry ( GlobalOrdinal  globalBlockRow,
GlobalOrdinal  globalBlockCol,
LocalOrdinal  blkRowSize,
LocalOrdinal  blkColSize,
LocalOrdinal  LDA,
const Teuchos::ArrayView< const Scalar > &  blockEntry 
) [inline]

Copy the contents of the input block-entry into the matrix.

This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.

If the specified block-entry already exists in the matrix, it will be over-written (replaced) by the input block-entry.

This method may be called any time (before or after fillComplete()).

Definition at line 601 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::sumIntoGlobalBlockEntry ( GlobalOrdinal  globalBlockRow,
GlobalOrdinal  globalBlockCol,
LocalOrdinal  blkRowSize,
LocalOrdinal  blkColSize,
LocalOrdinal  LDA,
const Teuchos::ArrayView< const Scalar > &  blockEntry 
) [inline]

Add the contents of the input block-entry into the matrix.

This method will create the specified block-entry if it doesn't already exist, but only if fillComplete() has not yet been called.

If the specified block-entry already exists in the matrix, the contents of the input block-entry will be added to the values that are already present.

This method may be called any time (before or after fillComplete()).

Definition at line 618 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getGlobalBlockEntryView ( GlobalOrdinal  globalBlockRow,
GlobalOrdinal  globalBlockCol,
LocalOrdinal &  numPtRows,
LocalOrdinal &  numPtCols,
Teuchos::ArrayRCP< const Scalar > &  blockEntry 
) const [inline]

Returns a const read-only view of a block-entry.

The arguments numPtRows and numPtCols are set to the dimensions of the block- entry on output. The stride (LDA in Blas terminology) is equal to numPtRows.

This method may be called any time (before or after fillComplete()), but will throw an exception if the specified block-entry doesn't already exist.

Definition at line 379 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getGlobalBlockEntryViewNonConst ( GlobalOrdinal  globalBlockRow,
GlobalOrdinal  globalBlockCol,
LocalOrdinal &  numPtRows,
LocalOrdinal &  numPtCols,
Teuchos::ArrayRCP< Scalar > &  blockEntry 
) [inline]

Returns a non-const read-write view of a block-entry.

Creates the block-entry if it doesn't already exist, and if:

  • the arguments numPtRows and numPtCols are set on entry (and nonzero),
  • and if fillComplete() has not yet been called.

Important Note: Be very careful managing the lifetime of this view. If fillComplete() has been called, and if you are running on a GPU, this view may be a copy of memory from the GPU, and your changes to the view won't be copied back to the GPU until your ArrayRCP is destroyed or set to Teuchos::null.

Definition at line 314 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getLocalBlockEntryView ( LocalOrdinal  localBlockRow,
LocalOrdinal  localBlockCol,
LocalOrdinal &  numPtRows,
LocalOrdinal &  numPtCols,
Teuchos::ArrayRCP< const Scalar > &  blockEntry 
) const [inline]

Returns a const read-only view of a block-entry.

The arguments numPtRows and numPtCols are set to the dimensions of the block- entry on output. The stride (LDA in Blas terminology) is equal to numPtRows. Throws an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't exist.

This method may only be called after fillComplete() has been called, and will throw an exception if the specified block-entry doesn't already exist.

Definition at line 472 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::getLocalBlockEntryViewNonConst ( LocalOrdinal  localBlockRow,
LocalOrdinal  localBlockCol,
LocalOrdinal &  numPtRows,
LocalOrdinal &  numPtCols,
Teuchos::ArrayRCP< Scalar > &  blockEntry 
) [inline]

Returns a non-const read-write view of a block-entry.

The arguments numPtRows and numPtCols are set to the dimensions of the block- entry on output. The stride (LDA in Blas terminology) is equal to numPtRows. Throws an exception if fillComplete() has not yet been called, or if the specified block-entry doesn't exist.

Important Note: Be very careful managing the lifetime of this view. If fillComplete() has been called, and if you are running on a GPU, this view may be a copy of memory from the GPU, and your changes to the view won't be copied back to the GPU until your ArrayRCP is destroyed or set to Teuchos::null.

This method may only be called after fillComplete() has been called, and will throw an exception if the specified block-entry doesn't already exist.

Definition at line 427 of file Tpetra_VbrMatrix_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node , class LocalMatOps >
void Tpetra::VbrMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node, LocalMatOps >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const [inline, virtual]

Print the object with some verbosity level to a FancyOStream object.

Reimplemented from Teuchos::Describable.

Definition at line 805 of file Tpetra_VbrMatrix_def.hpp.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Generated on Wed Apr 13 10:21:42 2011 for Tpetra Matrix/Vector Services by  doxygen 1.6.3