#include <Trilinos_Util_CrsMatrixGallery.h>
Inheritance diagram for Trilinos_Util::CrsMatrixGallery:
Extraction methods. | |
| Epetra_CrsMatrix * | GetMatrix () |
| Returns a pointer to the CrsMatrix. | |
| Epetra_CrsMatrix & | GetMatrixRef () |
| Epetra_Vector * | GetExactSolution () |
| Returns a pointer to the exact solution. | |
| Epetra_Vector * | GetStartingSolution () |
| Returns a pointer to the starting solution (typically, for HB problems). | |
| Epetra_Vector * | 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 (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 string &FileName, const bool UseSparse=true) |
| Print matrix on file in MATLAB format. | |
| ostream & | operator<< (ostream &os, const Trilinos_Util::CrsMatrixGallery &G) |
| Print out detailed information about the problem at hand. | |
Public Member Functions | |
Constructors/Destructor. | |
| CrsMatrixGallery (const string name, const Epetra_Comm &comm) | |
| Triutils_Gallery Constructor. | |
| CrsMatrixGallery (const string name, const Epetra_Map &map) | |
| Creates an Triutils_Gallery object using a given map. | |
| ~CrsMatrixGallery () | |
| Triutils_Gallery destructor. | |
Setting methods | |
| int | Set (const string parameter, const int value) |
| Sets a gallery options using an interger value. | |
| int | Set (const string parameter, const string value) |
| Sets a gallery options using a C++ string . | |
| int | Set (const string parameter, const double value) |
| Sets a gallery options using an double value. | |
| int | Set (const 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 | |
Creation methods. | |
| void | CreateMap () |
| Creates a map. | |
| void | CreateMatrix () |
| Creates the CrdMatrix. | |
| void | CreateExactSolution () |
| Creates the exact solution. | |
| void | CreateStartingSolution () |
| Creates the starting solution. | |
| void | CreateRHS () |
| Create the RHS corresponding to the desired exact solution. | |
| void | CreateEye () |
| void | CreateMatrixDiag () |
| void | CreateMatrixTriDiag () |
| void | CreateMatrixLaplace1d () |
| void | CreateMatrixLaplace1dNeumann () |
| void | CreateMatrixCrossStencil2d () |
| void | CreateMatrixCrossStencil2dVector () |
| void | CreateMatrixLaplace2d () |
| void | CreateMatrixLaplace2d_9pt () |
| void | CreateMatrixRecirc2d () |
| void | CreateMatrixRecirc2dDivFree () |
| void | CreateMatrixLaplace2dNeumann () |
| void | CreateMatrixUniFlow2d () |
| void | CreateMatrixLaplace3d () |
| void | CreateMatrixCrossStencil3d () |
| void | CreateMatrixCrossStencil3dVector () |
| void | CreateMatrixLehmer () |
| void | CreateMatrixMinij () |
| void | CreateMatrixRis () |
| void | CreateMatrixHilbert () |
| void | CreateMatrixJordblock () |
| void | CreateMatrixCauchy () |
| void | CreateMatrixFiedler () |
| void | CreateMatrixHanowa () |
| void | CreateMatrixKMS () |
| void | CreateMatrixParter () |
| void | CreateMatrixPei () |
| void | CreateMatrixOnes () |
| void | CreateMatrixVander () |
| void | ReadMatrix () |
| 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) |
| 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_Vector * | ExactSolution_ |
| Epetra_Vector * | StartingSolution_ |
| Epetra_Vector * | rhs_ |
| Epetra_Map * | map_ |
| Epetra_LinearProblem * | LinearProblem_ |
| string | name_ |
| int | NumGlobalElements_ |
| int | NumMyElements_ |
| int * | MyGlobalElements_ |
| string | MapType_ |
| string | ExactSolutionType_ |
| string | StartingSolutionType_ |
| string | ExpandType_ |
| string | RhsType_ |
| int | nx_ |
| int | ny_ |
| int | nz_ |
| int | mx_ |
| int | my_ |
| int | mz_ |
| int | NumPDEEqns_ |
| 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_ |
| string | FileName_ |
| string | ErrorMsg |
| string | OutputMsg |
| bool | verbose_ |
|
||||||||||||
|
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.
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...
|
|
||||||||||||
|
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.
|
|
|
Triutils_Gallery destructor.
|
|
|
Computes the 2-norm of the difference between the starting solution and the exact solution.
|
|
|
Computes the 2-norm of the residual.
|
|
|
Creates the exact solution.
|
|
|
|
|
|
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:
|
|
|
Creates the CrdMatrix.
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
Create the RHS corresponding to the desired exact solution.
|
|
|
Creates the starting solution.
|
|
||||||||||||||||||||||||||||||||
|
|
|
||||||||||||||||
|
|
|
||||||||||||||||
|
Get pointers to double vectors containing coordinates of points.
|
|
|
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:
|
|
|
Returns a pointer to Epetra_LinearProblem.
|
|
|
Returns a pointer the internally stored Map.
|
|
|
|
|
|
Returns a pointer to the CrsMatrix.
|
|
|
|
|
||||||||||||||||||||||||||||||||
|
|
|
||||||||||||||||||||||||||||||||||||||||||||
|
|
|
|
Returns a pointer to the rhs corresponding to the selected exact solution.
|
|
|
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
|
|
|
|
|
|
Print out matrix and vectors.
|
|
|
|
|
|
Sets gallery options using values passed from the shell.
|
|
||||||||||||
|
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 a gallery options using an double value.
|
|
||||||||||||
|
Sets a gallery options using a C++ string .
|
|
||||||||||||
|
Sets a gallery options using an interger value.
|
|
|
|
|
|
|
|
||||||||||||
|
Print matrix on file in MATLAB format.
|
|
|
|
|
||||||||||||
|
Print out detailed information about the problem at hand.
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
1.3.9.1