OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
MueLuBottomSolver.hpp
Go to the documentation of this file.
1//
2// Class MueLuBottomSolver
3// Interface to the SAAMG solver of MueLu.
4//
5// Copyright (c) 2017 - 2020, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
6// All rights reserved
7//
8// Implemented as part of the PhD thesis
9// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
10//
11// This file is part of OPAL.
12//
13// OPAL is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// You should have received a copy of the GNU General Public License
19// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
20//
21#include <AMReX.H>
23
24template <class Level>
26 const std::string& reuse)
27 : BottomSolver<Teuchos::RCP<amr::matrix_t>,
28 Teuchos::RCP<amr::multivector_t>,
29 Level>()
30 , A_mp(Teuchos::null)
31 , nSweeps_m(4)
32 , rebalance_m(rebalance)
33 , setupTimer_m(IpplTimings::getTimer("AMR MG bsolver setup"))
34{
35 this->initMueLuList_m(reuse);
36
37 factory_mp = Teuchos::rcp( new pListInterpreter_t(mueluList_m) );
38
39 // empty multigrid hierarchy with a finest level only
40 hierarchy_mp = factory_mp->CreateHierarchy();
41}
42
43
44template <class Level>
45void MueLuBottomSolver<Level>::solve(const Teuchos::RCP<mv_t>& x,
46 const Teuchos::RCP<mv_t>& b)
47{
48 // MueLu requires Xpetra multivectors (wrap them)
49 Teuchos::RCP<xmv_t> xx = MueLu::TpetraMultiVector_To_XpetraMultiVector(x);
50 Teuchos::RCP<xmv_t> xb = MueLu::TpetraMultiVector_To_XpetraMultiVector(b);
51
52 // InitialGuessIsZero = true
53 // startLevel = 0
54 hierarchy_mp->Iterate(*xb, *xx, nSweeps_m, true, 0);
55
56 // put multivector back
57 x->assign(*util_t::MV2NonConstTpetraMV2(*xx));
58}
59
60
61template <class Level>
62void MueLuBottomSolver<Level>::setOperator(const Teuchos::RCP<matrix_t>& A,
63 Level* level_p)
64{
65 IpplTimings::startTimer(setupTimer_m);
66
67 A_mp = MueLu::TpetraCrs_To_XpetraMatrix<scalar_t, lo_t, go_t, node_t>(A);
68 A_mp->SetFixedBlockSize(1); // only 1 DOF per node (pure Laplace problem)
69
70 Teuchos::RCP<mv_t> coords_p = Teuchos::rcp(
71 new amr::multivector_t(A->getDomainMap(), AMREX_SPACEDIM, false)
72 );
73
74 const scalar_t* domain = level_p->geom.ProbLo();
75 const scalar_t* dx = level_p->cellSize();
76 for (amrex::MFIter mfi(level_p->grids, level_p->dmap, true);
77 mfi.isValid(); ++mfi)
78 {
79 const AmrBox_t& tbx = mfi.tilebox();
80 const lo_t* lo = tbx.loVect();
81 const lo_t* hi = tbx.hiVect();
82
83 for (lo_t i = lo[0]; i <= hi[0]; ++i) {
84 for (lo_t j = lo[1]; j <= hi[1]; ++j) {
85#if AMREX_SPACEDIM == 3
86 for (lo_t k = lo[2]; k <= hi[2]; ++k) {
87#endif
88 AmrIntVect_t iv(D_DECL(i, j, k));
89 go_t gidx = level_p->serialize(iv);
90
91 coords_p->replaceGlobalValue(gidx, 0, domain[0] + (0.5 + i) * dx[0]);
92 coords_p->replaceGlobalValue(gidx, 1, domain[1] + (0.5 + j) * dx[1]);
93#if AMREX_SPACEDIM == 3
94 coords_p->replaceGlobalValue(gidx, 2, domain[2] + (0.5 + k) * dx[2]);
95 }
96#endif
97 }
98 }
99 }
100
101 Teuchos::RCP<xmv_t> coordinates = MueLu::TpetraMultiVector_To_XpetraMultiVector(coords_p);
102
103 Teuchos::RCP<mv_t> nullspace = Teuchos::rcp(new mv_t(A->getRowMap(), 1));
104 Teuchos::RCP<xmv_t> xnullspace = MueLu::TpetraMultiVector_To_XpetraMultiVector(nullspace);
105 xnullspace->putScalar(1.0);
106 hierarchy_mp->GetLevel(0)->Set("Nullspace", xnullspace);
107 hierarchy_mp->GetLevel(0)->Set("Coordinates", coordinates);
108 hierarchy_mp->setDefaultVerbLevel(Teuchos::VERB_HIGH);
109 hierarchy_mp->IsPreconditioner(false);
110 hierarchy_mp->GetLevel(0)->Set("A", A_mp);
111
112 Teuchos::RCP<level_t> finest_p = hierarchy_mp->GetLevel(0);
113 finest_p->Set("A", A_mp);
114
115 factory_mp->SetupHierarchy(*hierarchy_mp);
116
117 this->isInitialized_m = true;
118
119 IpplTimings::stopTimer(setupTimer_m);
120}
121
122
123template <class Level>
125 return nSweeps_m;
126}
127
128
129template <class Level>
130std::string
132
133 std::map<std::string, std::string> map;
134 map["NONE"] = "none";
135 map["RP"] = "RP";
136 map["RAP"] = "RAP";
137 map["SYMBOLIC"] = "S";
138 map["FULL"] = "full";
139
140 auto muelu = map.find(reuse);
141
142 if ( muelu == map.end() )
143 throw OpalException("MueLuBottomSolver::convertToMueLuReuseOption()",
144 "No MueLu reuse option '" + reuse + "'.");
145
146 return muelu->second;
147}
148
149
150template <class Level>
151void MueLuBottomSolver<Level>::initMueLuList_m(const std::string& reuse) {
152 mueluList_m.set("problem: type", "Poisson-3D");
153 mueluList_m.set("verbosity", "low");
154 mueluList_m.set("number of equations", 1);
155 mueluList_m.set("max levels", 8);
156 mueluList_m.set("cycle type", "V");
157
158 mueluList_m.set("coarse: max size", 200);
159 mueluList_m.set("multigrid algorithm", "sa");
160 mueluList_m.set("sa: damping factor", 1.33); // default: 1.33
161 mueluList_m.set("sa: use filtered matrix", true);
162 mueluList_m.set("filtered matrix: reuse eigenvalue", false); // false: more expensive
163
164 mueluList_m.set("repartition: enable", rebalance_m);
165 mueluList_m.set("repartition: rebalance P and R", rebalance_m);
166 mueluList_m.set("repartition: partitioner", "zoltan2");
167 mueluList_m.set("repartition: min rows per proc", 800);
168 mueluList_m.set("repartition: start level", 2);
169
170 Teuchos::ParameterList reparms;
171 reparms.set("algorithm", "multijagged"); // rcb
172 // reparms.set("partitioning_approach", "partition");
173
174 mueluList_m.set("repartition: params", reparms);
175
176 mueluList_m.set("smoother: type", "CHEBYSHEV");
177 mueluList_m.set("smoother: pre or post", "both");
178 Teuchos::ParameterList smparms;
179 smparms.set("chebyshev: degree", 3);
180 smparms.set("chebyshev: assume matrix does not change", false);
181 smparms.set("chebyshev: zero starting solution", true);
182 mueluList_m.set("smoother: params", smparms);
183
184 mueluList_m.set("coarse: type", "RELAXATION");
185 Teuchos::ParameterList cparms;
186 cparms.set("relaxation: type", "Gauss-Seidel");
187 cparms.set("relaxation: sweeps", nSweeps_m);
188 cparms.set("relaxation: zero starting solution", true);
189 cparms.set("relaxation: use l1", true);
190 cparms.set("relaxation: l1 eta", 1.5);
191 mueluList_m.set("coarse: params", cparms);
192
193// mueluList_m.set("coarse: type", "KLU2");
194
195 mueluList_m.set("aggregation: type", "uncoupled");
196 mueluList_m.set("aggregation: min agg size", 3);
197 mueluList_m.set("aggregation: max agg size", 27);
198
199 mueluList_m.set("transpose: use implicit", false);
200
201 mueluList_m.set("reuse: type", reuse);
202}
Some AMR types used a lot.
Definition: AmrDefs.h:33
Tpetra::MultiVector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > multivector_t
void setOperator(const Teuchos::RCP< matrix_t > &A, Level *level_p=nullptr)
amr::matrix_t matrix_t
MueLu::ParameterListInterpreter< scalar_t, lo_t, go_t, node_t > pListInterpreter_t
std::size_t getNumIters()
amr::global_ordinal_t go_t
amr::AmrBox_t AmrBox_t
amr::multivector_t mv_t
MueLuBottomSolver(const bool &rebalance, const std::string &reuse)
void solve(const Teuchos::RCP< mv_t > &x, const Teuchos::RCP< mv_t > &b)
amr::AmrIntVect_t AmrIntVect_t
Teuchos::RCP< manager_t > factory_mp
sets up hierarchy
static std::string convertToMueLuReuseOption(const std::string &reuse)
void initMueLuList_m(const std::string &reuse)
Teuchos::ParameterList mueluList_m
amr::local_ordinal_t lo_t
Teuchos::RCP< hierarchy_t > hierarchy_mp
manages the multigrid hierarchy
amr::scalar_t scalar_t
The base class for all OPAL exceptions.
Definition: OpalException.h:28
static void stopTimer(TimerRef t)
Definition: IpplTimings.h:192
static void startTimer(TimerRef t)
Definition: IpplTimings.h:187