#include <Trilinos_Util_CrsMatrixGallery.h>

Public Member Functions | |
| CrsMatrixGallery (const std::string name, const Epetra_Comm &comm, bool UseLongLong=false) | |
| Triutils_Gallery Constructor. | |
| CrsMatrixGallery (const std::string name, const Epetra_Map &map) | |
| Creates an Triutils_Gallery object using a given map. | |
| ~CrsMatrixGallery () | |
| Triutils_Gallery destructor. | |
| int | Set (const std::string parameter, const int value) |
| Sets a gallery options using an interger value. | |
| int | Set (const std::string parameter, const std::string value) |
| Sets a gallery options using a C++ string . | |
| int | Set (const std::string parameter, const double value) |
| Sets a gallery options using an double value. | |
| int | Set (const std::string parameter, const Epetra_Vector &value) |
| Sets a gallery options using an Epetra_Vector. | |
| int | Set (Trilinos_Util::CommandLineParser &CLP) |
| Sets gallery options using values passed from the shell. | |
Protected Member Functions | |
| template<typename int_type > | |
| int_type *& | MyGlobalElementsPtr () |
| template<typename int_type > | |
| std::vector< int_type > & | MapMapRef () |
| void | CreateMap () |
| Creates a map. | |
| template<typename int_type > | |
| void | TCreateMap () |
| void | CreateMatrix () |
| Creates the CrdMatrix. | |
| template<typename int_type > | |
| void | TCreateMatrix () |
| template<typename int_type > | |
| void | TCreateExactSolution () |
| Creates the exact solution. | |
| void | CreateExactSolution () |
| void | CreateStartingSolution () |
| Creates the starting solution. | |
| template<typename int_type > | |
| void | TCreateRHS () |
| Create the RHS corresponding to the desired exact solution. | |
| void | CreateRHS () |
| template<typename int_type > | |
| void | CreateEye () |
| template<typename int_type > | |
| void | CreateMatrixDiag () |
| template<typename int_type > | |
| void | CreateMatrixTriDiag () |
| template<typename int_type > | |
| void | CreateMatrixLaplace1d () |
| template<typename int_type > | |
| void | CreateMatrixLaplace1dNeumann () |
| template<typename int_type > | |
| void | CreateMatrixCrossStencil2d () |
| template<typename int_type > | |
| void | CreateMatrixCrossStencil2dVector () |
| template<typename int_type > | |
| void | CreateMatrixLaplace2d () |
| template<typename int_type > | |
| void | CreateMatrixLaplace2d_BC () |
| template<typename int_type > | |
| void | CreateMatrixLaplace2d_9pt () |
| template<typename int_type > | |
| void | CreateMatrixStretched2d () |
| template<typename int_type > | |
| void | CreateMatrixRecirc2d () |
| template<typename int_type > | |
| void | CreateMatrixRecirc2dDivFree () |
| template<typename int_type > | |
| void | CreateMatrixLaplace2dNeumann () |
| template<typename int_type > | |
| void | CreateMatrixUniFlow2d () |
| template<typename int_type > | |
| void | CreateMatrixLaplace3d () |
| template<typename int_type > | |
| void | CreateMatrixCrossStencil3d () |
| template<typename int_type > | |
| void | CreateMatrixCrossStencil3dVector () |
| template<typename int_type > | |
| void | CreateMatrixLehmer () |
| template<typename int_type > | |
| void | CreateMatrixMinij () |
| template<typename int_type > | |
| void | CreateMatrixRis () |
| template<typename int_type > | |
| void | CreateMatrixHilbert () |
| template<typename int_type > | |
| void | CreateMatrixJordblock () |
| template<typename int_type > | |
| void | CreateMatrixCauchy () |
| template<typename int_type > | |
| void | CreateMatrixFiedler () |
| template<typename int_type > | |
| void | CreateMatrixHanowa () |
| template<typename int_type > | |
| void | CreateMatrixKMS () |
| template<typename int_type > | |
| void | CreateMatrixParter () |
| template<typename int_type > | |
| void | CreateMatrixPei () |
| template<typename int_type > | |
| void | CreateMatrixOnes () |
| template<typename int_type > | |
| void | CreateMatrixVander () |
| template<typename int_type > | |
| void | TReadMatrix () |
| void | GetNeighboursCartesian2d (const int i, const int nx, const int ny, int &left, int &right, int &lower, int &upper) |
| void | GetNeighboursCartesian3d (const int i, const int nx, const int ny, const int nz, int &left, int &right, int &lower, int &upper, int &below, int &above) |
| template<typename int_type > | |
| void | TGetCartesianCoordinates (double *&x, double *&y, double *&z) |
| void | ZeroOutData () |
| void | SetupCartesianGrid2D () |
| void | SetupCartesianGrid3D () |
| void | ExactSolQuadXY (double x, double y, double &u) |
| void | ExactSolQuadXY (double x, double y, double &u, double &ux, double &uy, double &uxx, double &uyy) |
Protected Attributes | |
| const Epetra_Comm * | comm_ |
| Epetra_CrsMatrix * | matrix_ |
| Epetra_MultiVector * | ExactSolution_ |
| Epetra_MultiVector * | StartingSolution_ |
| Epetra_MultiVector * | rhs_ |
| Epetra_Map * | map_ |
| Epetra_LinearProblem * | LinearProblem_ |
| std::string | name_ |
| long long | NumGlobalElements_ |
| int | NumMyElements_ |
| int * | MyGlobalElements_int_ |
| std::vector< int > | MapMap_int_ |
| long long * | MyGlobalElements_LL_ |
| std::vector< long long > | MapMap_LL_ |
| std::string | MapType_ |
| bool | ContiguousMap_ |
| std::string | ExactSolutionType_ |
| std::string | StartingSolutionType_ |
| std::string | ExpandType_ |
| std::string | RhsType_ |
| int | nx_ |
| int | ny_ |
| int | nz_ |
| int | mx_ |
| int | my_ |
| int | mz_ |
| double | lx_ |
| double | ly_ |
| double | lz_ |
| int | NumPDEEqns_ |
| int | NumVectors_ |
| Epetra_Vector * | VectorA_ |
| Epetra_Vector * | VectorB_ |
| Epetra_Vector * | VectorC_ |
| Epetra_Vector * | VectorD_ |
| Epetra_Vector * | VectorE_ |
| Epetra_Vector * | VectorF_ |
| Epetra_Vector * | VectorG_ |
| double | a_ |
| double | b_ |
| double | c_ |
| double | d_ |
| double | e_ |
| double | f_ |
| double | g_ |
| double | alpha_ |
| double | beta_ |
| double | gamma_ |
| double | delta_ |
| double | conv_ |
| double | diff_ |
| double | source_ |
| double | epsilon_ |
| std::string | FileName_ |
| std::string | ErrorMsg |
| std::string | OutputMsg |
| bool | verbose_ |
| bool | UseLongLong_ |
| Epetra_CrsMatrix * | GetMatrix () |
| Returns a pointer to the CrsMatrix. | |
| Epetra_CrsMatrix & | GetMatrixRef () |
| Epetra_MultiVector * | GetExactSolution () |
| Returns a pointer to the exact solution. | |
| Epetra_MultiVector * | GetStartingSolution () |
| Returns a pointer to the starting solution (typically, for HB problems). | |
| Epetra_MultiVector * | GetRHS () |
| Returns a pointer to the rhs corresponding to the selected exact solution. | |
| const Epetra_Map * | GetMap () |
| Returns a pointer the internally stored Map. | |
| const Epetra_Map & | GetMapRef () |
| Epetra_LinearProblem * | GetLinearProblem () |
| Returns a pointer to Epetra_LinearProblem. | |
| void | ComputeResidual (double *residual) |
| Computes the 2-norm of the residual. | |
| void | ComputeDiffBetweenStartingAndExactSolutions (double *residual) |
| Computes the 2-norm of the difference between the starting solution and the exact solution. | |
| void | PrintMatrixAndVectors (std::ostream &os) |
| Print out matrix and vectors. | |
| void | PrintMatrixAndVectors () |
| void | GetCartesianCoordinates (double *&x, double *&y, double *&z) |
| Get pointers to double vectors containing coordinates of points. | |
| int | WriteMatrix (const std::string &FileName, const bool UseSparse=true) |
| Print matrix on file in MATLAB format. | |
| std::ostream & | operator<< (std::ostream &os, const Trilinos_Util::CrsMatrixGallery &G) |
| Print out detailed information about the problem at hand. | |
| Trilinos_Util::CrsMatrixGallery::CrsMatrixGallery | ( | const std::string | name, |
| const Epetra_Comm & | comm, | ||
| bool | UseLongLong = false |
||
| ) |
Triutils_Gallery Constructor.
Creates a Triutils_Gallery instance.
The first parameter is the name of the matrix. We refer to the Trilinos Tutorial for a detailed description of available matrices.
An example of program using this class is reported below.
int main(int argc, char *argv[]) { #ifdef HAVE_MPI MPI_Init(&argc,&argv); Epetra_MpiComm Comm (MPI_COMM_WORLD); #else Epetra_SerialComm Comm; #endif // create an Epetra matrix reading an H/B matrix Trilinos_Util_CrsMatrixGallery Gallery("hb", Comm); // set the name of the matrix Gallery.Set("matrix name", "bcsstk14.rsa"); Epetra_CrsMatrix * A; Epetra_Vector * ExactSolution; Epetra_Vector * RHS; Epetra_Vector * StartingSolution; // at this point the matrix is read from file A = Gallery.GetMatrix(); ExactSolution = Gallery.GetExactSolution(); // at this point the RHS is allocated and filled RHS = Gallery.GetRHS(); StartingSolution = Gallery.GetStartingSolution(); // create linear problem Epetra_LinearProblem Problem(A,StartingSolution,RHS); // create AztecOO instance AztecOO Solver(Problem); Solver.SetAztecOption( AZ_precond, AZ_dom_decomp ); Solver.Iterate(1000,1E-9); // compute residual double residual; Gallery.ComputeResidual(&residual); if( Comm.MyPID()==0 ) cout << "||b-Ax||_2 = " << residual << endl; Gallery.ComputeDiffBetweenStartingAndExactSolutions(&residual); if( Comm.MyPID()==0 ) cout << "||x_exact - x||_2 = " << residual << endl; #ifdef HAVE_MPI MPI_Finalize() ; #endif return 0 ; }
Class CommandLineParser can be used as well. In this case, one may decide to use the following:
Trilinos_Util::CommandLineParser CLP(argc,argv); // set a problem with no matrix name Trilinos_Util::CrsMatrixGallery Gallery("", Comm); // read parameters and settings from the shell line G.Set(CLP); // continue with your code...
| In | comm - Epetra communicator |
| Trilinos_Util::CrsMatrixGallery::CrsMatrixGallery | ( | const std::string | name, |
| const Epetra_Map & | map | ||
| ) |
Creates an Triutils_Gallery object using a given map.
Create a Triutils_Gallery object using an Epetra_Map. Problem size must match the elements in map.
| In | name - definition of the problem to be created. |
| In | map - Epetra_Map |
Triutils_Gallery destructor.
References ExactSolution_, LinearProblem_, map_, matrix_, rhs_, StartingSolution_, VectorA_, VectorB_, VectorC_, VectorD_, VectorE_, VectorF_, VectorG_, and ZeroOutData().
| void Trilinos_Util::CrsMatrixGallery::ComputeDiffBetweenStartingAndExactSolutions | ( | double * | residual | ) |
Computes the 2-norm of the difference between the starting solution and the exact solution.
| void Trilinos_Util::CrsMatrixGallery::ComputeResidual | ( | double * | residual | ) |
Computes the 2-norm of the residual.
| void Trilinos_Util::CrsMatrixGallery::CreateExactSolution | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateEye | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMap | ( | void | ) | [protected] |
Creates a map.
Creates an Epetra_Map. Before calling this function, the problem size must have been specified.
CreateMap() allows some different maps. The type of map is set using Set("map",value). Value is a string, defined as:
| void Trilinos_Util::CrsMatrixGallery::CreateMatrix | ( | void | ) | [protected] |
Creates the CrdMatrix.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixCauchy | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixCrossStencil2d | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixCrossStencil2dVector | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixCrossStencil3d | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixCrossStencil3dVector | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixDiag | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixFiedler | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixHanowa | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixHilbert | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixJordblock | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixKMS | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace1d | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace1dNeumann | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace2d | ( | void | ) | [protected] |
References Scaling.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace2d_9pt | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace2d_BC | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace2dNeumann | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLaplace3d | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixLehmer | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixMinij | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixOnes | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixParter | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixPei | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixRecirc2d | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixRecirc2dDivFree | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixRis | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixStretched2d | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixTriDiag | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixUniFlow2d | ( | void | ) | [protected] |
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::CreateMatrixVander | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateRHS | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::CreateStartingSolution | ( | void | ) | [protected] |
Creates the starting solution.
| void Trilinos_Util::CrsMatrixGallery::ExactSolQuadXY | ( | double | x, |
| double | y, | ||
| double & | u | ||
| ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::ExactSolQuadXY | ( | double | x, |
| double | y, | ||
| double & | u, | ||
| double & | ux, | ||
| double & | uy, | ||
| double & | uxx, | ||
| double & | uyy | ||
| ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::GetCartesianCoordinates | ( | double *& | x, |
| double *& | y, | ||
| double *& | z | ||
| ) |
Get pointers to double vectors containing coordinates of points.
| Epetra_MultiVector * Trilinos_Util::CrsMatrixGallery::GetExactSolution | ( | void | ) |
Returns a pointer to the exact solution.
Returns a pointer to the exact solution.
Some choices are available to define the exact solution, using Set("exact solution", value). value can be:
| Epetra_LinearProblem * Trilinos_Util::CrsMatrixGallery::GetLinearProblem | ( | void | ) |
Returns a pointer to Epetra_LinearProblem.
| const Epetra_Map * Trilinos_Util::CrsMatrixGallery::GetMap | ( | void | ) |
Returns a pointer the internally stored Map.
| const Epetra_Map & Trilinos_Util::CrsMatrixGallery::GetMapRef | ( | void | ) |
| Epetra_CrsMatrix * Trilinos_Util::CrsMatrixGallery::GetMatrix | ( | void | ) |
Returns a pointer to the CrsMatrix.
| Epetra_CrsMatrix & Trilinos_Util::CrsMatrixGallery::GetMatrixRef | ( | void | ) |
| void Trilinos_Util::CrsMatrixGallery::GetNeighboursCartesian2d | ( | const int | i, |
| const int | nx, | ||
| const int | ny, | ||
| int & | left, | ||
| int & | right, | ||
| int & | lower, | ||
| int & | upper | ||
| ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::GetNeighboursCartesian3d | ( | const int | i, |
| const int | nx, | ||
| const int | ny, | ||
| const int | nz, | ||
| int & | left, | ||
| int & | right, | ||
| int & | lower, | ||
| int & | upper, | ||
| int & | below, | ||
| int & | above | ||
| ) | [protected] |
| Epetra_MultiVector * Trilinos_Util::CrsMatrixGallery::GetRHS | ( | void | ) |
Returns a pointer to the rhs corresponding to the selected exact solution.
| Epetra_MultiVector * Trilinos_Util::CrsMatrixGallery::GetStartingSolution | ( | void | ) |
Returns a pointer to the starting solution (typically, for HB problems).
Returns a pointer to the starting solution. This is typically used while reading a HB problem. However, the user can set a starting solution using Set("starting solution", "value"). Value can be
| std::vector< int > & Trilinos_Util::CrsMatrixGallery::MapMapRef< int > | ( | ) | [inline, protected] |
| int *& Trilinos_Util::CrsMatrixGallery::MyGlobalElementsPtr< int > | ( | ) | [inline, protected] |
| void Trilinos_Util::CrsMatrixGallery::PrintMatrixAndVectors | ( | std::ostream & | os | ) |
Print out matrix and vectors.
| int Trilinos_Util::CrsMatrixGallery::Set | ( | const std::string | parameter, |
| const int | value | ||
| ) |
Sets a gallery options using an interger value.
| int Trilinos_Util::CrsMatrixGallery::Set | ( | const std::string | parameter, |
| const std::string | value | ||
| ) |
Sets a gallery options using a C++ string .
| int Trilinos_Util::CrsMatrixGallery::Set | ( | const std::string | parameter, |
| const double | value | ||
| ) |
Sets a gallery options using an double value.
| int Trilinos_Util::CrsMatrixGallery::Set | ( | const std::string | parameter, |
| const Epetra_Vector & | value | ||
| ) |
Sets a gallery options using an Epetra_Vector.
Sets a gallery options using an Epetra_Vector. The Epetra_Vector is copied into internal structures, and freed by the destructor.
Sets gallery options using values passed from the shell.
References Trilinos_Util_Map::Get(), and Trilinos_Util_Map::Has().
| void Trilinos_Util::CrsMatrixGallery::SetupCartesianGrid2D | ( | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::SetupCartesianGrid3D | ( | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::TCreateExactSolution | ( | void | ) | [protected] |
Creates the exact solution.
| void Trilinos_Util::CrsMatrixGallery::TCreateMap | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::TCreateMatrix | ( | void | ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::TCreateRHS | ( | void | ) | [protected] |
Create the RHS corresponding to the desired exact solution.
References UNDEF.
| void Trilinos_Util::CrsMatrixGallery::TGetCartesianCoordinates | ( | double *& | x, |
| double *& | y, | ||
| double *& | z | ||
| ) | [protected] |
| void Trilinos_Util::CrsMatrixGallery::TReadMatrix | ( | void | ) | [protected] |
| int Trilinos_Util::CrsMatrixGallery::WriteMatrix | ( | const std::string & | FileName, |
| const bool | UseSparse = true |
||
| ) |
Print matrix on file in MATLAB format.
| void Trilinos_Util::CrsMatrixGallery::ZeroOutData | ( | ) | [protected] |
References UNDEF.
Referenced by ~CrsMatrixGallery().
| std::ostream& operator<< | ( | std::ostream & | os, |
| const Trilinos_Util::CrsMatrixGallery & | G | ||
| ) | [friend] |
Print out detailed information about the problem at hand.
double Trilinos_Util::CrsMatrixGallery::a_ [protected] |
double Trilinos_Util::CrsMatrixGallery::alpha_ [protected] |
double Trilinos_Util::CrsMatrixGallery::b_ [protected] |
double Trilinos_Util::CrsMatrixGallery::beta_ [protected] |
double Trilinos_Util::CrsMatrixGallery::c_ [protected] |
const Epetra_Comm* Trilinos_Util::CrsMatrixGallery::comm_ [protected] |
Referenced by Trilinos_Util::operator<<().
bool Trilinos_Util::CrsMatrixGallery::ContiguousMap_ [protected] |
double Trilinos_Util::CrsMatrixGallery::conv_ [protected] |
double Trilinos_Util::CrsMatrixGallery::d_ [protected] |
double Trilinos_Util::CrsMatrixGallery::delta_ [protected] |
double Trilinos_Util::CrsMatrixGallery::diff_ [protected] |
double Trilinos_Util::CrsMatrixGallery::e_ [protected] |
double Trilinos_Util::CrsMatrixGallery::epsilon_ [protected] |
std::string Trilinos_Util::CrsMatrixGallery::ErrorMsg [protected] |
Epetra_MultiVector* Trilinos_Util::CrsMatrixGallery::ExactSolution_ [protected] |
Referenced by Trilinos_Util::operator<<(), and ~CrsMatrixGallery().
std::string Trilinos_Util::CrsMatrixGallery::ExactSolutionType_ [protected] |
Referenced by Trilinos_Util::operator<<().
std::string Trilinos_Util::CrsMatrixGallery::ExpandType_ [protected] |
double Trilinos_Util::CrsMatrixGallery::f_ [protected] |
std::string Trilinos_Util::CrsMatrixGallery::FileName_ [protected] |
double Trilinos_Util::CrsMatrixGallery::g_ [protected] |
double Trilinos_Util::CrsMatrixGallery::gamma_ [protected] |
Epetra_LinearProblem* Trilinos_Util::CrsMatrixGallery::LinearProblem_ [protected] |
Referenced by ~CrsMatrixGallery().
double Trilinos_Util::CrsMatrixGallery::lx_ [protected] |
double Trilinos_Util::CrsMatrixGallery::ly_ [protected] |
double Trilinos_Util::CrsMatrixGallery::lz_ [protected] |
Epetra_Map* Trilinos_Util::CrsMatrixGallery::map_ [protected] |
Referenced by ~CrsMatrixGallery().
std::vector<int> Trilinos_Util::CrsMatrixGallery::MapMap_int_ [protected] |
std::vector<long long> Trilinos_Util::CrsMatrixGallery::MapMap_LL_ [protected] |
std::string Trilinos_Util::CrsMatrixGallery::MapType_ [protected] |
Referenced by Trilinos_Util::operator<<().
Epetra_CrsMatrix* Trilinos_Util::CrsMatrixGallery::matrix_ [protected] |
Referenced by Trilinos_Util::operator<<(), and ~CrsMatrixGallery().
int Trilinos_Util::CrsMatrixGallery::mx_ [protected] |
int Trilinos_Util::CrsMatrixGallery::my_ [protected] |
int* Trilinos_Util::CrsMatrixGallery::MyGlobalElements_int_ [protected] |
long long* Trilinos_Util::CrsMatrixGallery::MyGlobalElements_LL_ [protected] |
int Trilinos_Util::CrsMatrixGallery::mz_ [protected] |
std::string Trilinos_Util::CrsMatrixGallery::name_ [protected] |
Referenced by Trilinos_Util::operator<<().
long long Trilinos_Util::CrsMatrixGallery::NumGlobalElements_ [protected] |
Referenced by Trilinos_Util::operator<<().
int Trilinos_Util::CrsMatrixGallery::NumMyElements_ [protected] |
int Trilinos_Util::CrsMatrixGallery::NumPDEEqns_ [protected] |
Referenced by Trilinos_Util::operator<<().
int Trilinos_Util::CrsMatrixGallery::NumVectors_ [protected] |
int Trilinos_Util::CrsMatrixGallery::nx_ [protected] |
int Trilinos_Util::CrsMatrixGallery::ny_ [protected] |
int Trilinos_Util::CrsMatrixGallery::nz_ [protected] |
std::string Trilinos_Util::CrsMatrixGallery::OutputMsg [protected] |
Epetra_MultiVector* Trilinos_Util::CrsMatrixGallery::rhs_ [protected] |
Referenced by Trilinos_Util::operator<<(), and ~CrsMatrixGallery().
std::string Trilinos_Util::CrsMatrixGallery::RhsType_ [protected] |
double Trilinos_Util::CrsMatrixGallery::source_ [protected] |
Epetra_MultiVector* Trilinos_Util::CrsMatrixGallery::StartingSolution_ [protected] |
Referenced by ~CrsMatrixGallery().
std::string Trilinos_Util::CrsMatrixGallery::StartingSolutionType_ [protected] |
bool Trilinos_Util::CrsMatrixGallery::UseLongLong_ [protected] |
Epetra_Vector* Trilinos_Util::CrsMatrixGallery::VectorA_ [protected] |
Referenced by ~CrsMatrixGallery().
Epetra_Vector * Trilinos_Util::CrsMatrixGallery::VectorB_ [protected] |
Referenced by ~CrsMatrixGallery().
Epetra_Vector * Trilinos_Util::CrsMatrixGallery::VectorC_ [protected] |
Referenced by ~CrsMatrixGallery().
Epetra_Vector * Trilinos_Util::CrsMatrixGallery::VectorD_ [protected] |
Referenced by ~CrsMatrixGallery().
Epetra_Vector * Trilinos_Util::CrsMatrixGallery::VectorE_ [protected] |
Referenced by ~CrsMatrixGallery().
Epetra_Vector * Trilinos_Util::CrsMatrixGallery::VectorF_ [protected] |
Referenced by ~CrsMatrixGallery().
Epetra_Vector * Trilinos_Util::CrsMatrixGallery::VectorG_ [protected] |
Referenced by ~CrsMatrixGallery().
bool Trilinos_Util::CrsMatrixGallery::verbose_ [protected] |
1.7.6.1