| 
    EpetraExt
    Development
    
   | 
 
EpetraExt_BlockDiagMatrix: A class for storing distributed block matrices. More...
#include <EpetraExt_BlockDiagMatrix.h>

Public Member Functions | |
| EpetraExt_BlockDiagMatrix (const Epetra_BlockMap &Map, bool zero_out=true) | |
| Constructor - This map is the map of the vector this can be applied to.  More... | |
| EpetraExt_BlockDiagMatrix (const EpetraExt_BlockDiagMatrix &Source) | |
| Copy constructor.  More... | |
| virtual | ~EpetraExt_BlockDiagMatrix () | 
| Destructor.  More... | |
| EpetraExt_BlockDiagMatrix & | operator= (const EpetraExt_BlockDiagMatrix &Source) | 
| = Operator.  More... | |
| double * | operator[] (int index) | 
| Block access function.  More... | |
| const double * | operator[] (int index) const | 
| Block access function.  More... | |
  Public Member Functions inherited from Epetra_BLAS | |
| Epetra_BLAS (void) | |
| Epetra_BLAS (const Epetra_BLAS &BLAS) | |
| virtual | ~Epetra_BLAS (void) | 
| float | ASUM (const int N, const float *X, const int INCX=1) const | 
| double | ASUM (const int N, const double *X, const int INCX=1) const | 
| float | DOT (const int N, const float *X, const float *Y, const int INCX=1, const int INCY=1) const | 
| double | DOT (const int N, const double *X, const double *Y, const int INCX=1, const int INCY=1) const | 
| float | NRM2 (const int N, const float *X, const int INCX=1) const | 
| double | NRM2 (const int N, const double *X, const int INCX=1) const | 
| void | SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const | 
| void | SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const | 
| void | COPY (const int N, const float *X, float *Y, const int INCX=1, const int INCY=1) const | 
| void | COPY (const int N, const double *X, double *Y, const int INCX=1, const int INCY=1) const | 
| int | IAMAX (const int N, const float *X, const int INCX=1) const | 
| int | IAMAX (const int N, const double *X, const int INCX=1) const | 
| void | AXPY (const int N, const float ALPHA, const float *X, float *Y, const int INCX=1, const int INCY=1) const | 
| void | AXPY (const int N, const double ALPHA, const double *X, double *Y, const int INCX=1, const int INCY=1) const | 
| void | GEMV (const char TRANS, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *X, const float BETA, float *Y, const int INCX=1, const int INCY=1) const | 
| void | GEMV (const char TRANS, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *X, const double BETA, double *Y, const int INCX=1, const int INCY=1) const | 
| void | GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const | 
| void | GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const | 
| void | SYMM (const char SIDE, const char UPLO, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const | 
| void | SYMM (const char SIDE, const char UPLO, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const | 
| void | TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const float ALPHA, const float *A, const int LDA, float *B, const int LDB) const | 
| void | TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const double ALPHA, const double *A, const int LDA, double *B, const int LDB) const | 
| void | SYRK (const char UPLO, const char TRANS, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float BETA, float *C, const int LDC) const | 
| void | SYRK (const char UPLO, const char TRANS, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double BETA, double *C, const int LDC) const | 
Attribute set methods | |
| virtual int | SetUseTranspose (bool) | 
| SetUseTranspose - not implemented.  More... | |
| virtual int | SetParameters (Teuchos::ParameterList &List) | 
| SetParameters.  More... | |
| virtual int | Compute () | 
| Computes the inverse / factorization if such is set on the list.  More... | |
Attribute access functions | |
| virtual const char * | Label () const | 
| Returns a character std::string describing the operator.  More... | |
| virtual bool | UseTranspose () const | 
| Returns the current UseTranspose setting.  More... | |
| virtual bool | HasNormInf () const | 
| Returns true if the this object can provide an approximate Inf-norm, false otherwise.  More... | |
| virtual const Epetra_Comm & | Comm () const | 
| Returns a pointer to the Epetra_Comm communicator associated with this operator.  More... | |
| virtual const Epetra_Map & | OperatorDomainMap () const | 
| Returns the Epetra_Map object associated with the domain of this operator.  More... | |
| virtual const Epetra_Map & | OperatorRangeMap () const | 
| Returns the Epetra_Map object associated with the range of this operator.  More... | |
| virtual const Epetra_BlockMap & | BlockMap () const | 
| Returns the Epetra_BlockMap object associated with the range of this operator.  More... | |
| double * | Values () const | 
| Returns a pointer to the array containing the blocks.  More... | |
| int | BlockSize (int LID) const | 
| Returns the size of the given block.  More... | |
| int | DataSize (int LID) const | 
| Returns the size of the data in the given block.  More... | |
| bool | ConstantBlockSize () const | 
| Returns true if the element size is constant.  More... | |
| int | NumMyBlocks () const | 
| Returns the number of local blocks.  More... | |
| int | NumGlobalBlocks () const | 
| Returns the number of global blocks.  More... | |
| long long | NumGlobalBlocks64 () const | 
| int | NumMyUnknowns () const | 
| Returns the number of local unknowns.  More... | |
| int | NumGlobalUnknowns () const | 
| Returns the number of global unknowns.  More... | |
| long long | NumGlobalUnknowns64 () const | 
| int | NumData () const | 
| Returns the size of the total Data block.  More... | |
| int | GetApplyMode () | 
| Gets apply mode info.  More... | |
| virtual void | Print (std::ostream &os) const | 
| Print method.  More... | |
Mathematical functions | |
| virtual int | Apply (const Epetra_MultiVector &, Epetra_MultiVector &) const | 
| Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.  More... | |
| virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const | 
| Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.  More... | |
| virtual double | NormInf () const | 
| NormInf - Not Implemented.  More... | |
| void | PutScalar (double value) | 
| PutScalar function.  More... | |
| virtual const Epetra_BlockMap & | DataMap () const | 
| Returns the Epetra_BlockMap object with the distribution of underlying values.  More... | |
Additional Inherited Members | 
EpetraExt_BlockDiagMatrix: A class for storing distributed block matrices.
A dense-block block-diagonal matrix with inversion/factorization capabilities.
This class has a rigid map structure — the Domain and Range maps must be the same.
Definition at line 66 of file EpetraExt_BlockDiagMatrix.h.
| EpetraExt_BlockDiagMatrix::EpetraExt_BlockDiagMatrix | ( | const Epetra_BlockMap & | Map, | 
| bool | zero_out = true  | 
        ||
| ) | 
Constructor - This map is the map of the vector this can be applied to.
Definition at line 56 of file EpetraExt_BlockDiagMatrix.cpp.
| EpetraExt_BlockDiagMatrix::EpetraExt_BlockDiagMatrix | ( | const EpetraExt_BlockDiagMatrix & | Source | ) | 
Copy constructor.
Definition at line 80 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  virtual | 
Destructor.
Definition at line 71 of file EpetraExt_BlockDiagMatrix.cpp.
| EpetraExt_BlockDiagMatrix & EpetraExt_BlockDiagMatrix::operator= | ( | const EpetraExt_BlockDiagMatrix & | Source | ) | 
= Operator.
| In | A - EpetraExt_BlockDiagMatrix to copy. | 
Definition at line 150 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  inline | 
Block access function.
Definition at line 92 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Block access function.
Definition at line 97 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
SetUseTranspose - not implemented.
Implements Epetra_Operator.
Definition at line 105 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  virtual | 
SetParameters.
Definition at line 127 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  virtual | 
Computes the inverse / factorization if such is set on the list.
Definition at line 175 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  inlinevirtual | 
Returns a character std::string describing the operator.
Implements Epetra_Operator.
Definition at line 120 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
Returns the current UseTranspose setting.
Implements Epetra_Operator.
Definition at line 123 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Implements Epetra_Operator.
Definition at line 126 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Implements Epetra_Operator.
Definition at line 129 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
Returns the Epetra_Map object associated with the domain of this operator.
Implements Epetra_Operator.
Definition at line 132 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
Returns the Epetra_Map object associated with the range of this operator.
Implements Epetra_Operator.
Definition at line 135 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inlinevirtual | 
Returns the Epetra_BlockMap object associated with the range of this operator.
Definition at line 138 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns a pointer to the array containing the blocks.
Definition at line 141 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the size of the given block.
Definition at line 144 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the size of the data in the given block.
Definition at line 147 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns true if the element size is constant.
Definition at line 150 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the number of local blocks.
Definition at line 153 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the number of global blocks.
Definition at line 157 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Definition at line 159 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the number of local unknowns.
Definition at line 162 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the number of global unknowns.
Definition at line 166 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Definition at line 168 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Returns the size of the total Data block.
Definition at line 171 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  inline | 
Gets apply mode info.
Definition at line 174 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  virtual | 
Print method.
Reimplemented from Epetra_DistObject.
Definition at line 307 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  inlinevirtual | 
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
| In | X - A Epetra_MultiVector of dimension NumVectors to multiply with matrix. | 
| Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. | 
Implements Epetra_Operator.
Definition at line 193 of file EpetraExt_BlockDiagMatrix.h.
      
  | 
  virtual | 
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.
| In | X - A Epetra_MultiVector of dimension NumVectors to solve for. | 
| Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. | 
Implements Epetra_Operator.
Definition at line 231 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  inlinevirtual | 
NormInf - Not Implemented.
Implements Epetra_Operator.
Definition at line 210 of file EpetraExt_BlockDiagMatrix.h.
| void EpetraExt_BlockDiagMatrix::PutScalar | ( | double | value | ) | 
PutScalar function.
Definition at line 143 of file EpetraExt_BlockDiagMatrix.cpp.
      
  | 
  inlinevirtual | 
Returns the Epetra_BlockMap object with the distribution of underlying values.
Definition at line 216 of file EpetraExt_BlockDiagMatrix.h.
 1.8.5