38 #ifndef _MLPRECONDITIONER_HPP_ 39 #define _MLPRECONDITIONER_HPP_ 41 #include <ml_MultiLevelPreconditioner.h> 43 #include <lifev/core/LifeV.hpp> 45 #include <lifev/core/filter/GetPot.hpp> 46 #include <lifev/core/array/MatrixEpetra.hpp> 47 #include <lifev/core/algorithm/Preconditioner.hpp> 48 #include <lifev/core/algorithm/PreconditionerIfpack.hpp> 59 class PreconditionerML:
67 typedef Preconditioner super;
69 typedef ML_Epetra::MultiLevelPreconditioner prec_raw_type;
70 typedef std::shared_ptr<prec_raw_type> prec_type;
72 typedef super::operator_raw_type operator_raw_type;
73 typedef super::operator_type operator_type;
82 PreconditionerML ( std::shared_ptr<Epetra_Comm> comm = std::shared_ptr<Epetra_Comm> (
new Epetra_MpiComm ( MPI_COMM_WORLD ) ) );
84 PreconditionerML ( std::shared_ptr<Epetra_Comm> comm = std::shared_ptr<Epetra_Comm> (
new Epetra_SerialComm ) );
88 virtual ~PreconditionerML();
94 PreconditionerML ( operator_type& matrix );
105 Int buildPreconditioner ( operator_type& matrix );
108 void resetPreconditioner();
117 virtual void createParametersList ( list_Type& list,
119 const std::string& section,
120 const std::string& subSection )
122 createMLList ( list, dataFile, section, subSection, M_comm->MyPID() == 0 );
132 static void createMLList ( list_Type& list,
134 const std::string& section,
135 const std::string& subSection =
"ML",
136 const bool& verbose =
true );
143 virtual Int ApplyInverse (
const Epetra_MultiVector& vector1, Epetra_MultiVector& vector2 )
const 145 return M_preconditioner->ApplyInverse ( vector1, vector2 );
153 virtual Int Apply (
const Epetra_MultiVector& vector1, Epetra_MultiVector& vector2 )
const 155 return M_preconditioner->Apply ( vector1, vector2 );
159 virtual void showMe ( std::ostream& output = std::cout )
const;
172 void setDataFromGetPot (
const GetPot& dataFile,
173 const std::string& section );
179 Int SetUseTranspose (
bool useTranspose =
false )
181 return M_preconditioner->SetUseTranspose (useTranspose);
201 void setVerticesCoordinates (std::shared_ptr<std::vector<Real> > xCoord,
202 std::shared_ptr<std::vector<Real> > yCoord,
203 std::shared_ptr<std::vector<Real> > zCoord);
215 super::prec_raw_type* preconditioner();
218 super::prec_type preconditionerPtr()
220 return M_preconditioner;
224 std::string preconditionerType()
232 return M_preconditioner->UseTranspose();
236 const Epetra_Map& OperatorRangeMap()
const 238 return M_preconditioner->OperatorRangeMap();
242 const Epetra_Map& OperatorDomainMap()
const 244 return M_preconditioner->OperatorDomainMap();
251 list_Type M_IfpackSubList;
252 std::shared_ptr<Epetra_Comm> M_comm;
256 operator_type M_operator;
258 prec_type M_preconditioner;
262 bool M_visualizationDataAvailable;
263 std::shared_ptr<std::vector<Real> > M_xCoord;
264 std::shared_ptr<std::vector<Real> > M_yCoord;
265 std::shared_ptr<std::vector<Real> > M_zCoord;
270 inline Preconditioner* createML()
272 return new PreconditionerML();
int32_type Int
Generic integer data.
static const LifeV::UInt elm_nodes_num[]
double Real
Generic real data.