1#include "Teko_MLPreconditionerFactory.hpp"
3#include "Teko_MLLinearOp.hpp"
6#include "ml_MultiLevelPreconditioner.h"
7#include "ml_epetra_utils.h"
8#include "ml_RowMatrix.h"
10#include "Thyra_get_Epetra_Operator.hpp"
18MLPreconditionerState::MLPreconditionerState() : isFilled_(false) {}
22 mlComm_ = Teuchos::rcpWithDealloc(comm,Teuchos::deallocFunctorDelete<ML_Comm>(cleanup_ML_Comm));
27 mlOp_ = Teuchos::rcpWithDealloc(op,Teuchos::deallocFunctorDelete<ML_Operator>(cleanup_ML_Operator));
40void MLPreconditionerState::cleanup_ML_Comm(ML_Comm * mlComm)
42 ML_Comm_Destroy(&mlComm);
45void MLPreconditionerState::cleanup_ML_Operator(ML_Operator * mlOp)
47 ML_Operator_Destroy(&mlOp);
55Teuchos::RCP<ML_Epetra::MultiLevelPreconditioner>
57 const std::vector<Teuchos::RCP<const Teuchos::ParameterList> > & coarseningParams)
60 std::vector<Teuchos::ParameterList> cpls(coarseningParams.size());
61 for(std::size_t i=0;i<coarseningParams.size();i++)
62 cpls[i] = *coarseningParams[i];
64 mlPreconditioner_ = rcp(
new ML_Epetra::MultiLevelPreconditioner(&*mlOp_, mainList, &diagonalOps_[0], &cpls[0], diagonalOps_.size()));
66 return mlPreconditioner_;
73MLPreconditionerFactory::MLPreconditionerFactory()
85 Teuchos::RCP<ML_Epetra::MultiLevelPreconditioner> mlPrecOp =
89 return Teuchos::rcp(
new MLLinearOp(mlPrecOp));
99 TEUCHOS_ASSERT(not mlState.
isFilled());
101 EpetraExt::CrsMatrix_SolverMap RowMatrixColMapTrans;
108 ML_Comm_Create(&mlComm);
111 ML_Operator * mlBlkMat = ML_Operator_Create(mlComm);
112 ML_Operator_BlkMatInit(mlBlkMat,mlComm,rowCnt,colCnt,ML_DESTROY_EVERYTHING);
114 std::vector<Epetra_RowMatrix *> aggMats;
115 for(
int r=0;r<rowCnt;r++) {
116 for(
int c=0;c<colCnt;c++) {
117 ML_Operator * tmp = 0;
118 Epetra_RowMatrix * EMat = 0;
119 std::stringstream ss;
121 ss <<
"fine_"<< r <<
"_" << c;
124 Teuchos::RCP<Epetra_CrsMatrix> crsMat
125 = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(Thyra::get_Epetra_Operator(*blo->getNonconstBlock(r,c)));
126 EMat = ML_Epetra::ModifyEpetraMatrixColMap(*crsMat,RowMatrixColMapTrans, ss.str().c_str());
128 aggMats.push_back(EMat);
131 tmp = ML_Operator_Create(mlComm);
132 ML_Operator_WrapEpetraMatrix(EMat, tmp);
135 ML_Operator_BlkMatInsert(mlBlkMat, tmp, r,c);
138 ML_Operator_BlkMatFinalize(mlBlkMat);
157 Teko_DEBUG_SCOPE(
"MLPreconditionerFactory::initializeFromParameterList",10);
160 coarseningParams_.clear();
163 blockRowCount_ = settings.get<
int>(
"Block Row Count");
167 mainParams_ = settings.sublist(
"Smoothing Parameters");
168 mainParams_.set<Teuchos::RCP<const Teko::InverseLibrary> >(
"smoother: teko inverse library",
getInverseLibrary());
172 const Teuchos::ParameterList & aggSublist = settings.sublist(
"Block Aggregation");
174 for(
int block=0;block<blockRowCount_;block++) {
176 std::stringstream ss;
177 ss <<
"Block " << block;
178 std::string sublistName = ss.str();
181 RCP<Teuchos::ParameterList> userSpec = rcp(
new Teuchos::ParameterList(aggSublist.sublist(sublistName)));
182 coarseningParams_.push_back(userSpec);
int blockRowCount(const BlockedLinearOp &blo)
Get the row count in a block linear operator.
An implementation of a state object for block preconditioners.
void initializeFromParameterList(const Teuchos::ParameterList &settings)
This function builds the internals of the preconditioner factory from a parameter list.
void fillMLPreconditionerState(const BlockedLinearOp &blo, MLPreconditionerState &mlState) const
Fills an ML preconditioner state object.
virtual Teuchos::RCP< PreconditionerState > buildPreconditionerState() const
Function that permits the construction of an arbitrary PreconditionerState object.
virtual LinearOp buildPreconditionerOperator(BlockedLinearOp &blo, BlockPreconditionerState &state) const
Function that is called to build the preconditioner for the linear operator that is passed in.
Contains operator internals need for ML.
bool isFilled() const
Has this object been filled yet.
void setAggregationMatrices(const std::vector< Epetra_RowMatrix * > &diags)
Set matrices to build multigrid hierarcy from.
void setIsFilled(bool value)
Set if ML internals been constructed yet.
void setMLOperator(ML_Operator *op)
set ML Operator pointer...it will be automatically cleaned up
void setMLComm(ML_Comm *comm)
set ML Comm pointer...it will be automatically cleaned up
Teuchos::RCP< ML_Epetra::MultiLevelPreconditioner > constructMLPreconditioner(const Teuchos::ParameterList &mainList, const std::vector< Teuchos::RCP< const Teuchos::ParameterList > > &coarseningParams)
Build an ML preconditioner object using the set of coarsening parameters.
Teuchos::RCP< const InverseLibrary > getInverseLibrary() const
Get the inverse library used by this preconditioner factory.