MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_EminPFactory_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_EMINPFACTORY_DEF_HPP
11#define MUELU_EMINPFACTORY_DEF_HPP
12
13#include <Xpetra_Matrix.hpp>
14#include <Xpetra_StridedMapFactory.hpp>
15
17#include "MueLu_ProductOperator.hpp"
18#include "MueLu_FlatOperator.hpp"
19
20#include "MueLu_Constraint.hpp"
21#include "MueLu_MasterList.hpp"
22#include "MueLu_Monitor.hpp"
23#include "MueLu_PerfUtils.hpp"
24#include "MueLu_Behavior.hpp"
25
26#ifdef HAVE_MUELU_BELOS
27#include <BelosLinearProblem.hpp>
28#include <BelosSolverFactory.hpp>
29#include <BelosXpetraAdapter.hpp>
30#endif
31
32namespace MueLu {
33
34template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35typename Teuchos::ScalarTraits<Scalar>::magnitudeType
37 ComputeProlongatorEnergyNorm(RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& A,
38 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& P,
39 Teuchos::FancyOStream& out) {
40 using MagnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
41 RCP<Matrix> AP;
42 AP = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A, false, *P, false, AP, out, true, true);
43 RCP<Matrix> PtAP;
44 PtAP = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*P, true, *AP, false, PtAP, out, true, true);
45 RCP<Vector> diag = Xpetra::VectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(PtAP->getRowMap());
46 PtAP->getLocalDiagCopy(*diag);
47 MagnitudeType norm = diag->norm1();
48 norm = Teuchos::ScalarTraits<MagnitudeType>::squareroot(norm);
49 return norm;
50}
51
52template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
54 RCP<ParameterList> validParamList = rcp(new ParameterList());
55
56#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
57 SET_VALID_ENTRY("emin: num iterations");
58 SET_VALID_ENTRY("emin: num reuse iterations");
59 SET_VALID_ENTRY("emin: iterative method");
60 {
61 validParamList->getEntry("emin: iterative method").setValidator(rcp(new Teuchos::StringValidator(Teuchos::tuple<std::string>("cg", "sd", "gmres"))));
62 }
63#undef SET_VALID_ENTRY
64
65 validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory for the matrix A used during internal iterations");
66 validParamList->set<RCP<const FactoryBase>>("P", Teuchos::null, "Generating factory for the initial guess");
67 validParamList->set<RCP<const FactoryBase>>("Constraint", Teuchos::null, "Generating factory for constraints");
68
69 validParamList->set<RCP<Matrix>>("P0", Teuchos::null, "Initial guess at P");
70 validParamList->set<bool>("Keep P0", false, "Keep an initial P0 (for reuse)");
71
72 validParamList->set<RCP<Constraint>>("Constraint0", Teuchos::null, "Initial Constraint");
73 validParamList->set<bool>("Keep Constraint0", false, "Keep an initial Constraint (for reuse)");
74
75 return validParamList;
76}
77
78template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80 Input(fineLevel, "A");
81
82 static bool isAvailableP0 = false;
83 static bool isAvailableConstraint0 = false;
84
85 // Here is a tricky little piece of code
86 // We don't want to request (aka call Input) when we reuse and P0 is available
87 // However, we cannot run something simple like this:
88 // if (!coarseLevel.IsAvailable("P0", this))
89 // Input(coarseLevel, "P");
90 // The reason is that it works fine during the request stage, but fails
91 // in the release stage as we _construct_ P0 during Build process. Therefore,
92 // we need to understand whether we are in Request or Release mode
93 // NOTE: This is a very unique situation, please try not to propagate the
94 // mode check any further
95
96 // #define externalSuppliedP0
97 if (coarseLevel.GetRequestMode() == Level::REQUEST) {
98#ifdef externalSuppliedP0
99 isAvailableP0 = coarseLevel.IsAvailable("P0");
100#else
101 isAvailableP0 = coarseLevel.IsAvailable("P0", this);
102#endif
103 isAvailableConstraint0 = coarseLevel.IsAvailable("Constraint0", this);
104 }
105
106 if (isAvailableP0 == false)
107 Input(coarseLevel, "P");
108
109 if (isAvailableConstraint0 == false)
110 Input(coarseLevel, "Constraint");
111}
112
113template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115 BuildP(fineLevel, coarseLevel);
116}
117
118template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
120 FactoryMonitor m(*this, "Prolongator minimization", coarseLevel);
121
122 const ParameterList& pL = GetParameterList();
123
124 // Get the matrix
125 RCP<Matrix> A = Get<RCP<Matrix>>(fineLevel, "A");
126
127 if (restrictionMode_) {
128 SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
129
130 A = Utilities::Transpose(*A, true);
131 }
132
133 // Get/make constraint operator
134 RCP<Constraint> X;
135 if (coarseLevel.IsAvailable("Constraint0", this)) {
136 // Reuse data
137 X = coarseLevel.Get<RCP<Constraint>>("Constraint0", this);
138 GetOStream(Runtime0) << "Reusing Constraint0" << std::endl;
139
140 } else {
141 // Construct data
142 X = Get<RCP<Constraint>>(coarseLevel, "Constraint");
143 }
144
145 // Get/make initial guess
146 // NOTE: the main assumption here that P0 satisfies both constraints:
147 // - nonzero pattern
148 // - nullspace preservation
149 RCP<Matrix> P0;
150 int numIts;
151#ifdef externalSuppliedP0
152 if (coarseLevel.IsAvailable("P0")) {
153#else
154 if (coarseLevel.IsAvailable("P0", this)) {
155#endif
156 // Reuse data
157#ifdef externalSuppliedP0
158 P0 = coarseLevel.Get<RCP<Matrix>>("P0");
159#else
160 P0 = coarseLevel.Get<RCP<Matrix>>("P0", this);
161#endif
162 numIts = pL.get<int>("emin: num reuse iterations");
163 GetOStream(Runtime0) << "Reusing P0" << std::endl;
164
165 } else {
166 GetOStream(Runtime0) << "Getting P from coarseLevel" << std::endl;
167 // Construct data
168 P0 = Get<RCP<Matrix>>(coarseLevel, "P");
169 numIts = pL.get<int>("emin: num iterations");
170 }
171
172 if (Behavior::debug() && IsPrint(Statistics1)) {
173 SubFactoryMonitor m2(*this, "Statistics", coarseLevel);
174 GetOStream(Statistics1) << "Energy norm of P0: " << ComputeProlongatorEnergyNorm(A, P0, GetOStream(Statistics0)) << std::endl;
175 GetOStream(Statistics1) << "Constraint residual norm of P0: " << X->ResidualNorm(P0) << std::endl;
176 }
177
178 std::string solverType = pL.get<std::string>("emin: iterative method");
179 if (solverType == "cg")
180 solverType = "Pseudo Block CG";
181 else if (solverType == "gmres")
182 solverType = "Pseudo Block GMRES";
183
184 RCP<Matrix> P;
185#ifdef HAVE_MUELU_BELOS
186 if (numIts > 0) {
187 // Construct diagonal preconditioner
188 RCP<const Vector> invDiagonal = Utilities::GetMatrixDiagonalInverse(*A);
189 auto invD = MatrixFactory::Build(invDiagonal);
190 auto flatInvD = rcp(new FlatOperator(invD, X));
191
192 // Construct projection of A
193 // In principle we should be using the operator proj * A * proj.
194 // But since the initial guess already satisfies the constraints and all subsequent iterates do too,
195 // the application of proj on the right side is a no-op.
196 auto flatA = rcp(new FlatOperator(A, X));
197 std::vector<RCP<Operator>> ops = {X, flatA}; // Instead of ops = {X, flatA, X};
198 RCP<Operator> projectedFlatA = rcp(new ProductOperator(ops));
199
200 // Vectors for rhs and solution
201 auto B = MultiVectorFactory::Build(projectedFlatA->getRangeMap(), 1);
202 auto vecP = MultiVectorFactory::Build(projectedFlatA->getDomainMap(), 1);
203 // initialize solution vector to P0
204 X->AssignMatrixEntriesToVector(*P0, *vecP);
205
206 auto belosProjectedFlatA = rcp(new Belos::XpetraOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(projectedFlatA));
207 auto belosFlatInvDiag = rcp(new Belos::XpetraOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(flatInvD));
208
209 {
210 using MV = Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
211 using OP = Belos::OperatorT<MV>;
212
213 auto problem = rcp(new Belos::LinearProblem<Scalar, MV, OP>());
214 problem->setOperator(belosProjectedFlatA);
215 problem->setRightPrec(belosFlatInvDiag);
216 problem->setLHS(vecP);
217 problem->setRHS(B);
218 problem->setLabel("Emin");
219 TEUCHOS_ASSERT(problem->setProblem());
220
221 auto belosList = rcp(new Teuchos::ParameterList());
222 belosList->set("Timer Label", "Emin");
223 belosList->set("Maximum Iterations", numIts);
224 belosList->set("Implicit Residual Scaling", "None");
225 belosList->set("Verbosity", Belos::Errors + Belos::Warnings + Belos::StatusTestDetails);
226 belosList->set("Output Frequency", 1);
227 belosList->set("Output Style", Belos::Brief);
228 auto out = GetMueLuOStream();
229 belosList->set("Output Stream", out->getOStream());
230
231 Belos::SolverFactory<Scalar, MV, OP> solverFactory;
232 auto solver = solverFactory.create(solverType, belosList);
233 solver->setProblem(problem);
234
235 if (numIts > 0) solver->solve();
236 }
237 // Convert from vector to matrix
238 P = X->GetMatrixWithEntriesFromVector(*vecP);
239 } else {
240 P = Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(P0);
241 }
242 if (P0->IsView("stridedMaps"))
243 P->CreateView("stridedMaps", P0);
244#else
245 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Energy minimization multigrid requires Belos to be enabled.");
246#endif
247
248 if (Behavior::debug() && IsPrint(Statistics1)) {
249 SubFactoryMonitor m2(*this, "Statistics", coarseLevel);
250 // Compute constraint residual norm if we are not in a debug build
251 GetOStream(Statistics1) << "Energy norm of P: " << ComputeProlongatorEnergyNorm(A, P, GetOStream(Statistics0)) << std::endl;
252 GetOStream(Statistics1) << "Norm of constraint residual: " << X->ResidualNorm(P) << std::endl;
253 }
254
255 // NOTE: EXPERIMENTAL and FRAGILE
256 if (!P->IsView("stridedMaps")) {
257 if (A->IsView("stridedMaps") == true) {
258 GetOStream(Runtime1) << "Using A to fillComplete P" << std::endl;
259
260 // FIXME: X->GetPattern() actually returns a CrsGraph.
261 // CrsGraph has no knowledge of Xpetra's sup/Matrix views. As such,
262 // it has no idea about strided maps. We create one, which is
263 // most likely incorrect for many use cases.
264 std::vector<size_t> stridingInfo(1, 1);
265 RCP<const StridedMap> dMap = StridedMapFactory::Build(X->GetPattern()->getDomainMap(), stridingInfo);
266
267 P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), dMap);
268
269 } else {
270 P->CreateView("stridedMaps", P->getRangeMap(), P->getDomainMap());
271 }
272 }
273
274 // Level Set
275 if (!restrictionMode_) {
276 // The factory is in prolongation mode
277 Set(coarseLevel, "P", P);
278
279 if (pL.get<bool>("Keep P0")) {
280 // NOTE: we must do Keep _before_ set as the Needs class only sets if
281 // a) data has been requested (which is not the case here), or
282 // b) data has some keep flag
283 coarseLevel.Keep("P0", this);
284 Set(coarseLevel, "P0", P);
285 }
286 if (pL.get<bool>("Keep Constraint0")) {
287 // NOTE: we must do Keep _before_ set as the Needs class only sets if
288 // a) data has been requested (which is not the case here), or
289 // b) data has some keep flag
290 coarseLevel.Keep("Constraint0", this);
291 Set(coarseLevel, "Constraint0", X);
292 }
293
294 if (IsPrint(Statistics2)) {
295 RCP<ParameterList> params = rcp(new ParameterList());
296 params->set("printLoadBalancingInfo", true);
297 params->set("printCommInfo", true);
298 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*P, "P", params);
299 }
300
301 } else {
302 // The factory is in restriction mode
303 RCP<Matrix> R;
304 {
305 SubFactoryMonitor m2(*this, "Transpose P", coarseLevel);
306
307 R = Utilities::Transpose(*P, true);
308 }
309
310 Set(coarseLevel, "R", R);
311
312 if (IsPrint(Statistics2)) {
313 RCP<ParameterList> params = rcp(new ParameterList());
314 params->set("printLoadBalancingInfo", true);
315 params->set("printCommInfo", true);
316 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*R, "R", params);
317 }
318 }
319}
320
321} // namespace MueLu
322
323#endif // MUELU_EMINPFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
static bool debug()
Whether MueLu is in debug mode.
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
static Teuchos::ScalarTraits< Scalar >::magnitudeType ComputeProlongatorEnergyNorm(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &P, Teuchos::FancyOStream &out)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Interprets a matrix as an operator that acts on a vector of nonzeros via SpGEMM.
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
RequestMode GetRequestMode() const
void Keep(const std::string &ename, const FactoryBase *factory)
Request to keep variable 'ename' generated by 'factory' after the setup phase.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
Takes a sequence of operators and applies their product.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool doLumped=false)
Extract Matrix Diagonal.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Statistics1
Print more statistics.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
@ Statistics0
Print statistics that do not involve significant additional computation.