Teko Version of the Day
Loading...
Searching...
No Matches
Teko_MLPreconditionerFactory.cpp
1// @HEADER
2// *****************************************************************************
3// Teko: A package for block and physics based preconditioning
4//
5// Copyright 2010 NTESS and the Teko contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#include "Teko_MLPreconditionerFactory.hpp"
11
12#include "Teko_MLLinearOp.hpp"
13
14#include "ml_include.h"
15#include "ml_MultiLevelPreconditioner.h"
16#include "ml_epetra_utils.h"
17#include "ml_RowMatrix.h"
18
19#include "Thyra_get_Epetra_Operator.hpp"
20
21namespace Teko {
22
24// MLPreconditionerState
26
27MLPreconditionerState::MLPreconditionerState() : isFilled_(false) {}
28
29void MLPreconditionerState::setMLComm(ML_Comm *comm) {
30 mlComm_ = Teuchos::rcpWithDealloc(comm, Teuchos::deallocFunctorDelete<ML_Comm>(cleanup_ML_Comm));
31}
32
33void MLPreconditionerState::setMLOperator(ML_Operator *op) {
34 mlOp_ =
35 Teuchos::rcpWithDealloc(op, Teuchos::deallocFunctorDelete<ML_Operator>(cleanup_ML_Operator));
36}
37
38void MLPreconditionerState::setIsFilled(bool value) { isFilled_ = value; }
39
40bool MLPreconditionerState::isFilled() const { return isFilled_; }
41
42void MLPreconditionerState::cleanup_ML_Comm(ML_Comm *mlComm) { ML_Comm_Destroy(&mlComm); }
43
44void MLPreconditionerState::cleanup_ML_Operator(ML_Operator *mlOp) { ML_Operator_Destroy(&mlOp); }
45
46void MLPreconditionerState::setAggregationMatrices(const std::vector<Epetra_RowMatrix *> &diags) {
47 diagonalOps_ = diags;
48}
49
50Teuchos::RCP<ML_Epetra::MultiLevelPreconditioner> MLPreconditionerState::constructMLPreconditioner(
51 const Teuchos::ParameterList &mainList,
52 const std::vector<Teuchos::RCP<const Teuchos::ParameterList> > &coarseningParams) {
53 TEUCHOS_ASSERT(isFilled());
54 std::vector<Teuchos::ParameterList> cpls(coarseningParams.size());
55 for (std::size_t i = 0; i < coarseningParams.size(); i++) cpls[i] = *coarseningParams[i];
56
57 mlPreconditioner_ = rcp(new ML_Epetra::MultiLevelPreconditioner(
58 &*mlOp_, mainList, &diagonalOps_[0], &cpls[0], diagonalOps_.size()));
59
60 return mlPreconditioner_;
61}
62
64// MLPreconditionerFactory
66
67MLPreconditionerFactory::MLPreconditionerFactory() {}
68
69LinearOp MLPreconditionerFactory::buildPreconditionerOperator(
70 BlockedLinearOp &blo, BlockPreconditionerState &state) const {
71 MLPreconditionerState &mlState = Teuchos::dyn_cast<MLPreconditionerState>(state);
72
73 if (not mlState.isFilled()) fillMLPreconditionerState(blo, mlState);
74 TEUCHOS_ASSERT(mlState.isFilled());
75
76 Teuchos::RCP<ML_Epetra::MultiLevelPreconditioner> mlPrecOp =
77 mlState.constructMLPreconditioner(mainParams_, coarseningParams_);
78
79 // return Thyra::epetraLinearOp(mlPrecOp);
80 return Teuchos::rcp(new MLLinearOp(mlPrecOp));
81}
82
83Teuchos::RCP<PreconditionerState> MLPreconditionerFactory::buildPreconditionerState() const {
84 return Teuchos::rcp(new MLPreconditionerState());
85}
86
87void MLPreconditionerFactory::fillMLPreconditionerState(const BlockedLinearOp &blo,
88 MLPreconditionerState &mlState) const {
89 TEUCHOS_ASSERT(not mlState.isFilled());
90
91 EpetraExt::CrsMatrix_SolverMap RowMatrixColMapTrans;
92
93 int rowCnt = blockRowCount(blo);
94 int colCnt = blockRowCount(blo);
95
96 // construct comm object...add to state
97 ML_Comm *mlComm;
98 ML_Comm_Create(&mlComm);
99 mlState.setMLComm(mlComm);
100
101 ML_Operator *mlBlkMat = ML_Operator_Create(mlComm);
102 ML_Operator_BlkMatInit(mlBlkMat, mlComm, rowCnt, colCnt, ML_DESTROY_EVERYTHING);
103
104 std::vector<Epetra_RowMatrix *> aggMats;
105 for (int r = 0; r < rowCnt; r++) {
106 for (int c = 0; c < colCnt; c++) {
107 ML_Operator *tmp = 0;
108 Epetra_RowMatrix *EMat = 0;
109 std::stringstream ss;
110
111 ss << "fine_" << r << "_" << c;
112
113 // extract crs matrix
114 Teuchos::RCP<Epetra_CrsMatrix> crsMat = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(
115 Thyra::get_Epetra_Operator(*blo->getNonconstBlock(r, c)));
116 EMat = ML_Epetra::ModifyEpetraMatrixColMap(*crsMat, RowMatrixColMapTrans, ss.str().c_str());
117 if (r == c) // setup diagonal scaling matrices
118 aggMats.push_back(EMat);
119
120 // extract ml sub operator
121 tmp = ML_Operator_Create(mlComm);
122 ML_Operator_WrapEpetraMatrix(EMat, tmp);
123
124 // add it to the block
125 ML_Operator_BlkMatInsert(mlBlkMat, tmp, r, c);
126 }
127 }
128 ML_Operator_BlkMatFinalize(mlBlkMat);
129
130 // finish setting up state object
131 mlState.setMLOperator(mlBlkMat);
132
133 // now set aggregation matrices
134 mlState.setAggregationMatrices(aggMats);
135
136 mlState.setIsFilled(true); // register state object as filled
137}
138
142void MLPreconditionerFactory::initializeFromParameterList(const Teuchos::ParameterList &settings) {
143 using Teuchos::RCP;
144 using Teuchos::rcp;
145
146 Teko_DEBUG_SCOPE("MLPreconditionerFactory::initializeFromParameterList", 10);
147
148 // clear initial state
149 coarseningParams_.clear();
150 blockRowCount_ = 0;
151
152 blockRowCount_ = settings.get<int>("Block Row Count");
153
154 // read in main parameter list: with smoothing information
156 mainParams_ = settings.sublist("Smoothing Parameters");
157 mainParams_.set<Teuchos::RCP<const Teko::InverseLibrary> >("smoother: teko inverse library",
158 getInverseLibrary());
159
160 // read in aggregation sub lists
162 const Teuchos::ParameterList &aggSublist = settings.sublist("Block Aggregation");
163
164 for (int block = 0; block < blockRowCount_; block++) {
165 // write out sub list name: "Block #"
166 std::stringstream ss;
167 ss << "Block " << block;
168 std::string sublistName = ss.str();
169
170 // grab sublist
171 RCP<Teuchos::ParameterList> userSpec =
172 rcp(new Teuchos::ParameterList(aggSublist.sublist(sublistName)));
173 coarseningParams_.push_back(userSpec);
174 }
175}
176
177} // end namespace Teko
int blockRowCount(const BlockedLinearOp &blo)
Get the row count in a block linear operator.
An implementation of a state object for block preconditioners.
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.