26 #ifndef OOMPH_TRILINOS_OPERATORS_HEADER
27 #define OOMPH_TRILINOS_OPERATORS_HEADER
31 #include <oomph-lib-config.h>
42 #include "ml_include.h"
43 #include "ml_MultiLevelPreconditioner.h"
111 void setup(Epetra_CrsMatrix* epetra_matrix_pt);
134 Epetra_CrsMatrix* epetra_matrix_pt) = 0;
220 ML_parameters.set(
"smoother: damping factor", smoother_damping);
238 ML_parameters.set(
"smoother: type",
"symmetric Gauss-Seidel");
A vector in the mathematical sense, initially developed for linear algebra type applications....
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
An interface to the Trilinos IFPACK class- provides a function to construct an IFPACK object,...
double & ilut_fill_level()
Access function for ILUT_fill_level.
virtual ~TrilinosIFPACKPreconditioner()
Destructor – empty, cleanup is done in base class.
double Absolute_threshold
Value of absolute threshold, used to peturb diagonal.
int & ilu_fill_level()
Access function for ILU_fill_level.
double & absolute_threshold()
Access function for the absolute threshold.
TrilinosIFPACKPreconditioner()
Constructor.
string Preconditioner_type
Type of ILU preconditioner.
int ILU_fill_level
Level of fill for "ILU".
void set_preconditioner_ILU()
Broken assignment operator.
void setup_trilinos_preconditioner(Epetra_CrsMatrix *epetra_matrix_pt)
Function to set up an IFPACK preconditioner. It is assumed Trilinos_matrix_pt points to a suitable ma...
void set_preconditioner_ILUT()
Function to set Preconditioner_type to "ILUT".
double ILUT_fill_level
Level of fill for "ILUT".
double Relative_threshold
Value of relative threshold, used to pertub diagonal.
int Overlap
Value of overlap level - used in parallel ILU.
double & relative_threshold()
Access function for the relative threshold.
TrilinosIFPACKPreconditioner(const TrilinosIFPACKPreconditioner &)=delete
Broken copy constructor.
An interface to the Trilinos ML class - provides a function to construct a serial ML object,...
void set_smoother_sweeps(int smoother_sweeps)
Function to set Smoother_sweeps.
void set_smoother_gauss_seidel()
Function to set smoother type to "symmetric Gauss-Seidel".
void set_max_levels(int max_levels)
Function to set maximum number of levels.
virtual ~TrilinosMLPreconditioner()
Destructor empty – clean up is done in base class.
void set_n_cycles(int n_cycles)
Function to set the number of cycles used.
static int Default_n_cycles
Default number of V cycles (one to be consistent with previous default) (It's an int because Trilinos...
TrilinosMLPreconditioner()
Constructor. Build with Smooth Aggretation (SA) default settings, but our own default number of V cyc...
Teuchos::ParameterList ML_parameters
void set_NSSA_default_values()
Broken assignment operator.
void set_DDML_default_values()
Set control flags to values 3-level algebraic domain decomposition.
void set_smoother_damping(double smoother_damping)
Function to set Smoother_damping.
void set_output(int output)
Function to set output - controls level of information output by ML.
void set_SA_default_values()
Set control flags to values for classical smoothed aggregation preconditioning.
void set_DD_default_values()
Set control flags to values for classical smoothed aggregation- based 2-level domain decomposition.
TrilinosMLPreconditioner(const TrilinosMLPreconditioner &)=delete
Broken copy constructor.
void setup_trilinos_preconditioner(Epetra_CrsMatrix *epetra_matrix_pt)
Function to set up the ML preconditioner. It is assumed Trilinos_matrix_pt points to a suitable matri...
void set_smoother_jacobi()
Function to set smoother type to "Jacobi".
Base class for Trilinos preconditioners as oomph-lib preconditioner.
Epetra_Operator * Epetra_preconditioner_pt
The preconditioner which will be set up using function setup_trilinos_preconditioner(....
Epetra_Operator * epetra_operator_pt() const
Access function to Epetra_preconditioner_pt (const version) For use with TrilinosAztecOOSolver.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
applies the preconditioner
Epetra_Operator *& epetra_operator_pt()
Access function to Epetra_preconditioner_pt. For use with TrilinosAztecOOSolver.
Epetra_CrsMatrix * Epetra_matrix_pt
Pointer used to store the epetra matrix - only used when this preconditioner is setup using the oomph...
virtual ~TrilinosPreconditionerBase()
Destructor.
void clean_up_memory()
deletes the preconditioner, matrices and maps
static double Cumulative_preconditioner_solve_time
Static double that accumulates the preconditioner solve time of all instantiations of this class....
void setup()
Broken assignment operator.
virtual void setup_trilinos_preconditioner(Epetra_CrsMatrix *epetra_matrix_pt)=0
Function to set up a specific Trilinos preconditioner. This is called by setup(......
TrilinosPreconditionerBase()
Constructor.
TrilinosPreconditionerBase(const TrilinosPreconditionerBase &)=delete
Broken copy constructor.
void output()
Doc the command line arguments.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...