MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MultiPhys_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_MULTIPHYS_DEF_HPP
11#define MUELU_MULTIPHYS_DEF_HPP
12
13#include <sstream>
14#include "MueLu_ConfigDefs.hpp"
15
16#include "Xpetra_Map.hpp"
18#include "Xpetra_MatrixUtils.hpp"
19
21#include "MueLu_SaPFactory.hpp"
22#include "MueLu_AggregationExportFactory.hpp"
23#include "MueLu_Utilities.hpp"
24#include "MueLu_Level.hpp"
25#include "MueLu_Hierarchy.hpp"
26#include "MueLu_RAPFactory.hpp"
27#include "MueLu_Monitor.hpp"
28#include "MueLu_PerfUtils.hpp"
29#include "MueLu_ParameterListInterpreter.hpp"
30#include "MueLu_HierarchyManager.hpp"
31#include <MueLu_HierarchyUtils.hpp>
35
36#ifdef HAVE_MUELU_CUDA
37#include "cuda_profiler_api.h"
38#endif
39
40namespace MueLu {
41
42template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> MultiPhys<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getDomainMap() const {
44 return AmatMultiphysics_->getDomainMap();
45}
46
47template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
48const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> MultiPhys<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getRangeMap() const {
49 return AmatMultiphysics_->getRangeMap();
50}
51
52template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
54 // this operator only makes sense for the Combo when using TransP for R
55
56 list.set("multigrid algorithm", "combine");
57 list.set("combine: numBlks", nBlks_);
58
59 // Make sure verbosity gets passed to the sublists
60 std::string verbosity = list.get("verbosity", "high");
62
63 arrayOfParamLists_.resize(nBlks_);
64 for (int i = 0; i < nBlks_; i++) {
65 std::string listName = "subblockList" + Teuchos::toString(i);
66 if (list.isSublist(listName)) {
67 arrayOfParamLists_[i] = Teuchos::rcpFromRef(list.sublist(listName));
68 } else
69 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Must provide sublist " + listName);
70
71 arrayOfParamLists_[i]->set("verbosity", arrayOfParamLists_[i]->get("verbosity", verbosity));
72 if (OmitSubblockSmoother_) {
73 arrayOfParamLists_[i]->set("smoother: pre or post", "none");
74 arrayOfParamLists_[i]->set("smoother: type", "none");
75 }
76 }
77
78 // Are we using Kokkos?
79 useKokkos_ = !Node::is_serial;
80 useKokkos_ = list.get("use kokkos refactor", useKokkos_);
81
82 paramListMultiphysics_ = Teuchos::rcpFromRef(list);
83}
84
85template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87 /*
88
89 Create a set of AMG hierarchies whose interpolation matrices are used to build on combined
90 AMG hierarchy for a multiphysics problem
91
92 */
93
94 //#ifdef HAVE_MUELU_CUDA
95 // if (paramListMultiphysics_.get<bool>("multiphysics: cuda profile setup", false)) cudaProfilerStart();
96 //#endif
97
98 std::string timerLabel;
99 if (reuse)
100 timerLabel = "MueLu MultiPhys: compute (reuse)";
101 else
102 timerLabel = "MueLu MultiPhys: compute";
103 RCP<Teuchos::TimeMonitor> tmCompute = getTimer(timerLabel);
104
106 // Generate the (iii,iii) Hierarchy
107
108 for (int iii = 0; iii < nBlks_; iii++) {
109 if (arrayOfCoords_ != Teuchos::null) {
110 if (arrayOfCoords_[iii] != Teuchos::null) {
111 arrayOfParamLists_[iii]->sublist("user data").set("Coordinates", arrayOfCoords_[iii]);
112 }
113 }
114
115 if (arrayOfMaterials_ != Teuchos::null) {
116 if (arrayOfMaterials_[iii] != Teuchos::null) {
117 arrayOfParamLists_[iii]->sublist("user data").set("Material", arrayOfMaterials_[iii]);
118 }
119 }
120
121 if (arrayOfNullspaces_ != Teuchos::null) {
122 if (arrayOfNullspaces_[iii] != Teuchos::null) {
123 arrayOfParamLists_[iii]->sublist("user data").set("Nullspace", arrayOfNullspaces_[iii]);
124 }
125 }
126
127 bool wantToRepartition = false;
128 if (paramListMultiphysics_->isParameter("repartition: enable"))
129 wantToRepartition = paramListMultiphysics_->get<bool>("repartition: enable");
130
131 arrayOfParamLists_[iii]->set("repartition: enable", wantToRepartition);
132 arrayOfParamLists_[iii]->set("repartition: rebalance P and R", true);
133 arrayOfParamLists_[iii]->set("repartition: explicit via new copy rebalance P and R", true);
134
135 if (paramListMultiphysics_->isParameter("repartition: use subcommunicators"))
136 arrayOfParamLists_[iii]->set("repartition: use subcommunicators", paramListMultiphysics_->isParameter("repartition: use subcommunicators"));
137 else
138 arrayOfParamLists_[iii]->set("repartition: use subcommunicators", true);
139 }
140 // repartitioning should only happen when createing the individual P's , not
141 // when combiing them
142
143 paramListMultiphysics_->set<bool>("repartition: enable", false);
144
145 bool useMaxLevels = false;
146 if (paramListMultiphysics_->isParameter("combine: useMaxLevels"))
147 useMaxLevels = paramListMultiphysics_->get<bool>("combine: useMaxLevels");
148
149 LO maxLevels = useMaxLevels ? 0 : std::numeric_limits<LO>::max();
150 for (int i = 0; i < nBlks_; i++) {
151 std::string operatorLabel = "MultiPhys (" + Teuchos::toString(i) + "," + Teuchos::toString(i) + ")";
152 arrayOfAuxMatrices_[i]->setObjectLabel(operatorLabel);
153 arrayOfHierarchies_[i] = MueLu::CreateXpetraPreconditioner(arrayOfAuxMatrices_[i], *arrayOfParamLists_[i]);
154 LO tempNlevels = arrayOfHierarchies_[i]->GetGlobalNumLevels();
155 if (useMaxLevels) {
156 if (tempNlevels > maxLevels) maxLevels = tempNlevels;
157 } else {
158 if (tempNlevels < maxLevels) maxLevels = tempNlevels;
159 }
160 }
161
162 hierarchyMultiphysics_ = rcp(new Hierarchy("Combo"));
163 for (LO i = 0; i < maxLevels; i++) {
164 hierarchyMultiphysics_->AddNewLevel();
165 }
166 for (int i = 0; i < nBlks_; i++) {
167 std::string subblkName = "Psubblock" + Teuchos::toString(i);
168 MueLu::HierarchyUtils<SC, LO, GO, NO>::CopyBetweenHierarchies(*(arrayOfHierarchies_[i]), *(hierarchyMultiphysics_), "P", subblkName, "RCP<Matrix>");
169
170 std::string subblkOpName = "Operatorsubblock" + Teuchos::toString(i);
171 MueLu::HierarchyUtils<SC, LO, GO, NO>::CopyBetweenHierarchies(*(arrayOfHierarchies_[i]), *(hierarchyMultiphysics_), "A", subblkOpName, "RCP<Matrix>");
172
173 // Copy remaining levels, if needed
174 if (useMaxLevels) {
175 const auto numLevelsBlk = arrayOfHierarchies_[i]->GetNumLevels();
176 const auto numGlobalLevelsBlk = arrayOfHierarchies_[i]->GetGlobalNumLevels();
177 if (numLevelsBlk == numGlobalLevelsBlk) {
178 auto crsLevel = arrayOfHierarchies_[i]->GetLevel(numLevelsBlk - 1);
179 TEUCHOS_ASSERT(crsLevel->IsAvailable("A"));
180 for (int levelId = numLevelsBlk; levelId < maxLevels; ++levelId) {
181 auto level = hierarchyMultiphysics_->GetLevel(levelId);
182 MueLu::HierarchyUtils<SC, LO, GO, NO>::CopyBetweenLevels(*crsLevel, *level, "A", subblkOpName, "RCP<Matrix>");
183 }
184 }
185 }
186 }
187 paramListMultiphysics_->set("coarse: max size", 1);
188 paramListMultiphysics_->set("max levels", maxLevels);
189
190 AmatMultiphysics_->setObjectLabel("A: block " + Teuchos::toString(nBlks_) + " x " + Teuchos::toString(nBlks_) + "multiphysics matrix");
191
192 // Rip off non-serializable data before validation
193 Teuchos::ParameterList nonSerialListMultiphysics, processedListMultiphysics;
194 MueLu::ExtractNonSerializableData(*paramListMultiphysics_, processedListMultiphysics, nonSerialListMultiphysics);
195
196 // Rip off the subblock List stuff as we don't need it any more and I think it messes up validator
197
198 Teuchos::ParameterList stripped;
199 for (ParameterList::ConstIterator inListEntry = processedListMultiphysics.begin(); inListEntry != processedListMultiphysics.end(); inListEntry++) {
200 const std::string& levelName = inListEntry->first;
201 if (levelName.find("subblockList") != 0) stripped.setEntry(inListEntry->first, inListEntry->second);
202 }
203
204 RCP<HierarchyManager<SC, LO, GO, NO>> mueLuFactory = rcp(new ParameterListInterpreter<SC, LO, GO, NO>(stripped, AmatMultiphysics_->getDomainMap()->getComm()));
205 hierarchyMultiphysics_->setlib(AmatMultiphysics_->getDomainMap()->lib());
206 hierarchyMultiphysics_->SetProcRankVerbose(AmatMultiphysics_->getDomainMap()->getComm()->getRank());
207
208 // We don't need nullspace or coordinates, since we don't use them when just combining prolongators that have been already created
209 hierarchyMultiphysics_->GetLevel(0)->Set("A", AmatMultiphysics_);
210
211 // Stick the non-serializible data on the hierarchy.
212 // Not sure that we need this, since we don't use it in building the multiphysics hierarchy
213 HierarchyUtils<SC, LO, GO, NO>::AddNonSerializableDataToHierarchy(*mueLuFactory, *hierarchyMultiphysics_, nonSerialListMultiphysics);
214 mueLuFactory->SetupHierarchy(*hierarchyMultiphysics_);
215
216 describe(GetOStream(Runtime0));
217}
218
219template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
220Teuchos::RCP<Teuchos::TimeMonitor> MultiPhys<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getTimer(std::string name, RCP<const Teuchos::Comm<int>> comm) const {
221 if (IsPrint(Timings)) {
222 if (!syncTimers_)
223 return Teuchos::rcp(new Teuchos::TimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name)));
224 else {
225 if (comm.is_null())
226 return Teuchos::rcp(new Teuchos::SyncTimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name), AmatMultiphysics_->getRowMap()->getComm().ptr()));
227 else
228 return Teuchos::rcp(new Teuchos::SyncTimeMonitor(*Teuchos::TimeMonitor::getNewTimer(name), comm.ptr()));
229 }
230 } else
231 return Teuchos::null;
232}
233
234template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
235void MultiPhys<Scalar, LocalOrdinal, GlobalOrdinal, Node>::resetMatrix(RCP<Matrix> AmatMultiphysics_new, bool ComputePrec) {
236 bool reuse = !AmatMultiphysics_.is_null();
237 AmatMultiphysics_ = AmatMultiphysics_new;
238 if (ComputePrec) compute(reuse);
239}
240
241template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
242void MultiPhys<Scalar, LocalOrdinal, GlobalOrdinal, Node>::applyInverse(const MultiVector& RHS, MultiVector& X) const {
243 hierarchyMultiphysics_->Iterate(RHS, X, 1, true);
244}
245
246template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
247void MultiPhys<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const MultiVector& RHS, MultiVector& X,
248 Teuchos::ETransp /* mode */,
249 Scalar /* alpha */,
250 Scalar /* beta */) const {
251 RCP<Teuchos::TimeMonitor> tm = getTimer("MueLu MultiPhys: solve");
252 hierarchyMultiphysics_->Iterate(RHS, X, 1, true);
253}
254
255template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
259
260template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
262 initialize(const Teuchos::RCP<Matrix>& AmatMultiPhysics,
263 const Teuchos::ArrayRCP<RCP<Matrix>> arrayOfAuxMatrices,
264 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfNullspaces,
265 const Teuchos::ArrayRCP<Teuchos::RCP<RealValuedMultiVector>> arrayOfCoords,
266 const int nBlks,
267 Teuchos::ParameterList& List,
268 const Teuchos::ArrayRCP<Teuchos::RCP<MultiVector>> arrayOfMaterials) {
269 arrayOfHierarchies_.resize(nBlks_);
270 for (int i = 0; i < nBlks_; i++) arrayOfHierarchies_[i] = Teuchos::null;
271
272 // Default settings
273 useKokkos_ = false;
274 enable_reuse_ = false;
275 syncTimers_ = false;
276
277 // set parameters
278 setParameters(List);
279
280} // initialize
281
282template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
284 describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel /* verbLevel */) const {
285 std::ostringstream oss;
286
287 RCP<const Teuchos::Comm<int>> comm = AmatMultiphysics_->getDomainMap()->getComm();
288
289 oss << "\n--------------------------------------------------------------------------------\n"
290 << "--- MultiPhysics Summary ---\n"
291 "--------------------------------------------------------------------------------"
292 << std::endl;
293 oss << std::endl;
294
295 GlobalOrdinal numRows;
296 GlobalOrdinal nnz;
297
298 AmatMultiphysics_->getRowMap()->getComm()->barrier();
299
300 for (int i = 0; i < nBlks_; i++) {
301 numRows = arrayOfAuxMatrices_[i]->getGlobalNumRows();
302 nnz = arrayOfAuxMatrices_[i]->getGlobalNumEntries();
303 Xpetra::global_size_t tt = numRows;
304 int rowspacer = 3;
305 while (tt != 0) {
306 tt /= 10;
307 rowspacer++;
308 }
309 tt = nnz;
310 int nnzspacer = 2;
311 while (tt != 0) {
312 tt /= 10;
313 nnzspacer++;
314 }
315
316 oss << "block " << std::setw(rowspacer) << " rows " << std::setw(nnzspacer) << " nnz " << std::setw(9) << " nnz/row" << std::endl;
317 oss << "(" << Teuchos::toString(i) << ", " << Teuchos::toString(i) << ")" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
318 }
319 oss << std::endl;
320
321 out << oss.str();
322
323 for (int i = 0; i < nBlks_; i++) {
324 arrayOfHierarchies_[i]->describe(out, GetVerbLevel());
325 }
326
327} // describe
328
329} // namespace MueLu
330
331#define MUELU_MULTIPHYS_SHORT
332#endif // ifdef MUELU_MULTIPHYS_DEF_HPP
Various adapters that will create a MueLu preconditioner that is an Xpetra::Matrix.
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception throws to report errors in the internal logical of the program.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
static void CopyBetweenLevels(Level &fromLevel, Level &toLevel, const std::string fromLabel, const std::string toLabel, const std::string dataType)
static void AddNonSerializableDataToHierarchy(HierarchyManager &HM, Hierarchy &H, const ParameterList &nonSerialList)
Add non-serializable data to Hierarchy.
static void CopyBetweenHierarchies(Hierarchy &fromHierarchy, Hierarchy &toHierarchy, const std::string fromLabel, const std::string toLabel, const std::string dataType)
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int > > comm=Teuchos::null) const
get a (synced) timer
void compute(bool reuse=false)
Setup the preconditioner.
void applyInverse(const MultiVector &RHS, MultiVector &X) const
apply standard MultiPhys cycle
void initialize(const Teuchos::RCP< Matrix > &AmatMultiPhysics, const Teuchos::ArrayRCP< RCP< Matrix > > arrayOfAuxMatrices, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfNullspaces, const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector > > arrayOfCoords, const int nBlks, Teuchos::ParameterList &List, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector > > arrayOfMaterials)
void setParameters(Teuchos::ParameterList &list)
Set parameters.
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
Namespace for MueLu classes and methods.
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Extract non-serializable data from level-specific sublists and move it to a separate parameter list.
@ Runtime0
One-liner description of what is happening.
@ Timings
Print all timing information.
MsgType toVerbLevel(const std::string &verbLevelStr)
Teuchos::RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > CreateXpetraPreconditioner(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > op, const Teuchos::ParameterList &inParamList)
Helper function to create a MueLu preconditioner that can be used by Xpetra.Given an Xpetra::Matrix,...