32 #ifndef TRILINOS_LINEAR_SOLVER_H 
   33 #define TRILINOS_LINEAR_SOLVER_H 
   51 #include "Epetra_MpiComm.h"  
   52 #include "Epetra_Map.h"  
   53 #include "Epetra_FEVbrMatrix.h"  
   54 #include "Epetra_FEVector.h" 
   55 #include "Epetra_FECrsGraph.h" 
   56 #include "Epetra_FECrsMatrix.h" 
   57 #include "Epetra_Import.h" 
   61 #include "AztecOO_StatusTestResNorm.h" 
   64 #include "ml_include.h" 
   65 #include "ml_MultiLevelPreconditioner.h" 
   69 #include "Ifpack_ConfigDefs.h" 
   70 #include "Ifpack_AdditiveSchwarz.h" 
   71 #include "Ifpack_ILUT.h" 
   72 #include "Ifpack_IC.h" 
   73 #include "Ifpack_ICT.h" 
   80 #define TRILINOS_CG_SOLVER 798 
   81 #define TRILINOS_GMRES_SOLVER 797 
   82 #define TRILINOS_BICGSTAB_SOLVER 795 
   85 #define NO_PRECONDITIONER 700 
   86 #define TRILINOS_DIAGONAL_PRECONDITIONER 702 
   87 #define TRILINOS_BLOCK_JACOBI_PRECONDITIONER 703 
   88 #define TRILINOS_ILU_PRECONDITIONER 704 
   89 #define TRILINOS_ILUT_PRECONDITIONER 705 
   90 #define TRILINOS_IC_PRECONDITIONER 706 
   91 #define TRILINOS_ICT_PRECONDITIONER 707 
   92 #define TRILINOS_ML_PRECONDITIONER 708 
   98   static Epetra_FEVector *
F;
 
   99   static Epetra_FEVbrMatrix *
K;
 
  100   static Epetra_Vector *
X;
 
  104   static Epetra_MpiComm *comm;
 
  105   static Epetra_FECrsGraph *K_graph;
 
  124   int Apply(
const Epetra_MultiVector &x, Epetra_MultiVector &y) 
const;
 
  130     return Trilinos::K->SetUseTranspose(use_transpose);
 
  134   int ApplyInverse(
const Epetra_MultiVector &X, Epetra_MultiVector &Y)
 const 
  164   const Epetra_Comm &
Comm()
 const 
  191           int& numGhostAndLocalNodes, 
int& nnz, 
const int *ltgSorted,
 
  192           const int *ltgUnsorted, 
const int *rowPtr, 
const int *colInd,
 
  193           int &
dof, 
int& cpp_index, 
int& proc_id);
 
  208           const double *lK, 
double *lR);
 
  211           double *x, 
const double *dirW, 
double &resNorm, 
double &initNorm,
 
  212           int &numIters, 
double &solverTime, 
double &dB, 
bool &converged,
 
  213           int &
lsType, 
double &relTol, 
int &maxIters, 
int &kspace,
 
  217           double &initNorm, 
int &numIters, 
double &solverTime,
 
  218           double &dB, 
bool &converged, 
int &
lsType, 
double &relTol,
 
  219           int &maxIters, 
int &kspace, 
int &precondType, 
bool &isFassem);
 
  228 void setPreconditioner(
int precondType, AztecOO &Solver);
 
  237               Epetra_Vector &diagonal);
 
  242 void printRHSToFile();
 
  244 void printSolutionToFile();
 
This class implements the pure virtual class Epetra_Operator for the AztecOO iterative solve which on...
Definition: trilinos_linear_solver.h:115
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Computes A_inv*x.
Definition: trilinos_linear_solver.h:134
int Apply(const Epetra_MultiVector &x, Epetra_MultiVector &y) const
Definition: trilinos_linear_solver.cpp:112
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with teh range of this operator.
Definition: trilinos_linear_solver.h:176
const char * Label() const
Returns a character string describing the operator.
Definition: trilinos_linear_solver.h:146
bool HasNormInf() const
Returns true if this object can provide an approx Inf-norm false otherwise.
Definition: trilinos_linear_solver.h:158
const Epetra_Comm & Comm() const
Returns pointer to Epetra_Comm communicator associated with this operator.
Definition: trilinos_linear_solver.h:164
const Epetra_Map & OperatorDomainMap() const
Returns Epetra_Map object assoicated with domain of this operator.
Definition: trilinos_linear_solver.h:170
double NormInf() const
Infinity norm for global stiffness does not add in the boundary term.
Definition: trilinos_linear_solver.h:140
int SetUseTranspose(bool use_transpose)
Definition: trilinos_linear_solver.h:128
bool UseTranspose() const
Returns current UseTranspose setting.
Definition: trilinos_linear_solver.h:152
Linear system of equations solver type.
Definition: ComMod.h:621
Initialize all Epetra types we need separate from Fortran.
Definition: trilinos_linear_solver.h:96
static Epetra_Vector * ghostX
Solution vector with nodes owned by processor followed by its ghost nodes.
Definition: trilinos_linear_solver.h:101
static Epetra_Import * Importer
Import ghostMap into blockMap to create ghost map.
Definition: trilinos_linear_solver.h:102
static Epetra_FEVector * bdryVec
Contribution from coupled neumann boundary conditions.
Definition: trilinos_linear_solver.h:103
static Epetra_Vector * X
Solution vector consisting of unique nodes owned by the processor.
Definition: trilinos_linear_solver.h:100
static Epetra_BlockMap * blockMap
Unique block map consisting of nodes owned by each processor.
Definition: trilinos_linear_solver.h:97
static Epetra_FEVector * F
Global block force vector.
Definition: trilinos_linear_solver.h:98
static Epetra_FEVbrMatrix * K
Global block stiffness matrix.
Definition: trilinos_linear_solver.h:99
int dof
Nodal degrees of freedom.
Definition: trilinos_linear_solver.cpp:78
void trilinos_global_solve_(const double *Val, const double *RHS, double *x, const double *dirW, double &resNorm, double &initNorm, int &numIters, double &solverTime, double &dB, bool &converged, int &lsType, double &relTol, int &maxIters, int &kspace, int &precondType)
Definition: trilinos_linear_solver.cpp:443
void trilinos_lhs_free_()
Definition: trilinos_linear_solver.cpp:969
void constructJacobiScaling(const double *dirW, Epetra_Vector &diagonal)
Definition: trilinos_linear_solver.cpp:888
void trilinos_lhs_create_(int &numGlobalNodes, int &numLocalNodes, int &numGhostAndLocalNodes, int &nnz, const int *ltgSorted, const int *ltgUnsorted, const int *rowPtr, const int *colInd, int &dof, int &cpp_index, int &proc_id)
Give function definitions which will be called through fortran.
Definition: trilinos_linear_solver.cpp:155
void trilinos_doassem_(int &numNodesPerElement, const int *eqN, const double *lK, double *lR)
Definition: trilinos_linear_solver.cpp:342
void setIFPACKPrec(AztecOO &Solver)
Definition: trilinos_linear_solver.cpp:833
void trilinos_solve_(double *x, const double *dirW, double &resNorm, double &initNorm, int &numIters, double &solverTime, double &dB, bool &converged, int &lsType, double &relTol, int &maxIters, int &kspace, int &precondType, bool &isFassem)
Definition: trilinos_linear_solver.cpp:546
void printMatrixToFile()
Definition: trilinos_linear_solver.cpp:1015
void trilinos_bc_create_(const double *v, bool &isCoupledBC)
Definition: trilinos_linear_solver.cpp:938
void setMLPrec(AztecOO &Solver)
Definition: trilinos_linear_solver.cpp:756
void checkDiagonalIsZero()
Definition: trilinos_linear_solver.cpp:864