26#ifndef OOMPH_TRILINOS_SOLVER_HEADER
27#define OOMPH_TRILINOS_SOLVER_HEADER
91 std::ostringstream error_message;
92 error_message <<
"SetUseTranspose() is a pure virtual Epetra_Operator "
93 <<
"member that is not required for a Preconditioner"
102 std::ostringstream error_message;
103 error_message <<
"Apply() is a pure virtual Epetra_Operator member"
104 <<
"that is not required for a Preconditioner" << std::endl;
132 unsigned nrow_local =
134 for (
unsigned i = 0;
i < nrow_local;
i++)
161 unsigned nrow_local =
163 for (
unsigned i = 0;
i < nrow_local;
i++)
176 std::ostringstream error_message;
177 error_message <<
"NormInf() is a pure virtual Epetra_Operator member"
178 <<
"that is not required for a Preconditioner" << std::endl;
192 std::ostringstream error_message;
194 <<
"UseTranspose() is a pure virtual Epetra_Operator member "
195 <<
"that is not required for a Preconditioner" << std::endl;
203 std::ostringstream error_message;
204 error_message <<
"HasNormInf() is a pure virtual Epetra_Operator member "
205 <<
"that is not required for a Preconditioner" << std::endl;
The conjugate gradient method.
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
Abstract base class for matrices of doubles – adds abstract interfaces for solving,...
A vector in the mathematical sense, initially developed for linear algebra type applications....
Base class for all linear iterative solvers. This merely defines standard interfaces for linear itera...
double Tolerance
Convergence tolerance.
unsigned Max_iter
Maximum number of iterations.
Preconditioner *& preconditioner_pt()
Access function to preconditioner.
OomphCommunicator * communicator_pt() const
const access to the communicator pointer
unsigned nrow_local() const
access function for the num of local rows on this processor. If no MPI then Nrow is returned.
bool Enable_resolve
Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be st...
An OomphLibError object which should be thrown when an run-time error is encountered....
An Epetra_Operator class for oomph-lib preconditioners. A helper class for TrilinosOomphLibPreconditi...
Epetra_SerialComm Operator_comm
An Epetra Serial Comm object.
Preconditioner * Oomph_lib_preconditioner_pt
A pointer to the oomph-lib preconditioner.
Epetra_MpiComm Operator_comm
An Epetra MPI Comm object.
OomphLibPreconditionerEpetraOperator(const OomphLibPreconditionerEpetraOperator &)=delete
Broken copy constructor.
bool Use_epetra_values
Use the epetra data within the vectors passed to the oomph-lib preconditioner. If this is true none o...
bool HasNormInf() const
Broken Epetra_Operator member - HasNormInf.
int ApplyInverse(const Epetra_MultiVector &epetra_r, Epetra_MultiVector &epetra_z) const
applies the oomph-lib preconditioner. Converts the Epetra vector applys the preconditioner by calling...
const char * Label() const
Epetra_Operator::Label - returns a string describing the operator.
void operator=(const OomphLibPreconditionerEpetraOperator &)=delete
Broken assignment operator.
std::string Preconditioner_label
a label for the preconditioner ( for Epetra_Operator::Label() )
bool UseTranspose() const
Broken Epetra_Operator member - UseTranspose.
const Epetra_Comm & Comm() const
Returns the Epetra MPI_Comm object.
double NormInf() const
Broken Epetra_Operator member - NormInf.
OomphLibPreconditionerEpetraOperator(Preconditioner *preconditioner_pt, bool use_epetra_values=false)
Constructor - takes the pointer to the oomph-lib preconditioner and the distribution of the precondit...
const Epetra_Map & OperatorDomainMap() const
Epetra_Operator member - OperatorDomainMap.
int SetUseTranspose(bool UseTranspose)
Broken Epetra_Operator member - SetUseTranspose.
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Broken Epetra_Operator member - Apply.
const Epetra_Map & OperatorRangeMap() const
Epetra_Operator member - OperatorRangeMap.
Epetra_Map * Operator_map_pt
A pointer to an Epetra_Map object - describes distribution of the preconditioner, in this instance it...
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
virtual void clean_up_memory()
Clean up memory (empty). Generic interface function.
virtual void preconditioner_solve(const DoubleVector &r, DoubleVector &z)=0
Apply the preconditioner. Pure virtual generic interface function. This method should apply the preco...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
An interface to the Trilinos AztecOO classes allowing it to be used as an Oomph-lib LinearSolver....
bool Using_problem_based_solve
Helper flag keeping track of whether we called the linear algebra or problem-based solve function.
void resolve(const DoubleVector &rhs, DoubleVector &solution)
Function to resolve a linear system using the existing solver data, allowing a solve with a new right...
bool If_oomphlib_preconditioner_use_epetra_values
if this solver is using an oomph-lib preconditioner then the vectors passed to preconditioner_solve(....
double & tolerance()
Access function to Tolerance.
void solver_setup(DoubleMatrixBase *const &matrix_pt)
Helper function for setting up the solver. Converts the oomph-lib matrices to Epetra matrices,...
void disable_aztecoo_workaround_for_epetra_matrix_setup()
Disable workaround for creating of epetra matrix that respects aztecoo's ordering requirements.
unsigned Iterations
Stores number of iterations used.
TrilinosAztecOOSolver()
Constructor.
Epetra_CrsMatrix * Epetra_matrix_pt
A pointer for the linear system matrix in Epetra_CrsMatrix format.
void clean_up_memory()
Clean up method - deletes the solver, the matrices and the preconditioner.
double jacobian_setup_time()
Function to return Jacobian_setup_time;.
bool Assemble_serial_jacobian
If true, when performing a problem based solve a serial matrix will be requested from Problem::get_ja...
void solve_using_AztecOO(Epetra_Vector *&rhs_pt, Epetra_Vector *&soln_pt)
Helper function performs the actual solve once the AztecOO solver is set up.
bool Delete_matrix
Trilinos copies matrix data from oomph-lib's own CRDoubleMatrix or DistributedCRDoubleMatrix to Trili...
void disable_delete_matrix()
Call if the matrix can not be deleted (default)
DoubleMatrixBase * Oomph_matrix_pt
Oomph lib matrix pointer.
Epetra_Operator * Epetra_preconditioner_pt
A pointer to the Epetra_Operator for the preconditioner. This is only used if the preconditioner NOT ...
~TrilinosAztecOOSolver()
Destructor - delete the solver and the matrices.
unsigned & solver_type()
Access function to Solver_type.
void disable_resolve()
Disable resolve function (overloads the LinearSolver disable_resolve function).
void solve(Problem *const &problem_pt, DoubleVector &solution)
Function which uses problem_pt's get_jacobian(...) function to generate a linear system which is then...
void enable_assemble_serial_jacobian()
Set the assembly of the serial jacobian when performing a problem-based solve.
double linear_solver_solution_time()
Function to return Linear_solver_solution_time.
AztecOO_solver_types
if this solver is using an oomph-lib preconditioner then the vectors passed to preconditioner_solve(....
unsigned & max_iter()
Access function to Max_iter.
void disable_assemble_serial_jacobian()
Unset the assembly of the serial jacobian.
TrilinosAztecOOSolver(const TrilinosAztecOOSolver &)=delete
Broken copy constructor.
void enable_aztecoo_workaround_for_epetra_matrix_setup()
Enable workaround for creating of epetra matrix that respects aztecoo's ordering requirements.
unsigned Solver_type
Defines which solver is set up - available types are defined in AztecOO_solver_types.
unsigned iterations() const
Acess function to Iterations.
double Linear_solver_solution_time
Stores time for the solution (excludes time to set up preconditioner)
Problem * Problem_pt
A pointer to the underlying problem (NULL if MATRIX based solve) The problem_pt is stored here in a p...
bool is_aztecoo_workaround_for_epetra_matrix_setup_enabled()
Is workaround for creating of epetra matrix that respects aztecoo's ordering requirements enabled?
void enable_delete_matrix()
Call if the matrix can be deleted.
AztecOO * AztecOO_solver_pt
Pointer to the AztecOO solver.
void operator=(const TrilinosAztecOOSolver &)=delete
Broken assignment operator.
bool Use_aztecoo_workaround_for_epetra_matrix_setup
Use workaround for creating of epetra matrix that respects aztecoo's ordering requirements.
double Jacobian_setup_time
Stores set up time for Jacobian.
Epetra_Map * create_epetra_map(const LinearAlgebraDistribution *const dist)
create an Epetra_Map corresponding to the LinearAlgebraDistribution
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).