MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_PgPFactory_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_PGPFACTORY_DEF_HPP
11#define MUELU_PGPFACTORY_DEF_HPP
12
13#include <vector>
14
15#include <Xpetra_Vector.hpp>
16#include <Xpetra_MultiVectorFactory.hpp>
17#include <Xpetra_VectorFactory.hpp>
18#include <Xpetra_Import.hpp>
19#include <Xpetra_ImportFactory.hpp>
20#include <Xpetra_Export.hpp>
21#include <Xpetra_ExportFactory.hpp>
22#include <Xpetra_Matrix.hpp>
23#include <Xpetra_MatrixMatrix.hpp>
24
26#include "MueLu_Monitor.hpp"
27#include "MueLu_PerfUtils.hpp"
29#include "MueLu_Utilities.hpp"
30
31namespace MueLu {
32
33template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35 RCP<ParameterList> validParamList = rcp(new ParameterList());
36
37 validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A used during the prolongator smoothing process");
38 validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Tentative prolongator factory");
39 validParamList->set<MinimizationNorm>("Minimization norm", DINVANORM, "Norm to be minimized");
40 validParamList->set<bool>("ReUseRowBasedOmegas", false, "Reuse omegas for prolongator for restrictor");
41
42 return validParamList;
43}
44
45template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
47 SetParameter("Minimization norm", ParameterEntry(minnorm)); // revalidate
48}
49
50template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
55
56template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
58 Input(fineLevel, "A");
59
60 // Get default tentative prolongator factory
61 // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
62 // -- Warning: Do not use directly initialPFact_. Use initialPFact instead everywhere!
63 RCP<const FactoryBase> initialPFact = GetFactory("P");
64 if (initialPFact == Teuchos::null) {
65 initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent");
66 }
67 coarseLevel.DeclareInput("P", initialPFact.get(), this);
68
69 /* If PgPFactory is reusing the row based damping parameters omega for
70 * restriction, it has to request the data here.
71 * we have the following scenarios:
72 * 1) Reuse omegas:
73 * PgPFactory.DeclareInput for prolongation mode requests A and P0
74 * PgPFactory.DeclareInput for restriction mode requests A, P0 and RowBasedOmega (call triggered by GenericRFactory)
75 * PgPFactory.Build for prolongation mode calculates RowBasedOmega and stores it (as requested)
76 * PgPFactory.Build for restriction mode reuses RowBasedOmega (and Releases the data with the Get call)
77 * 2) do not reuse omegas
78 * PgPFactory.DeclareInput for prolongation mode requests A and P0
79 * PgPFactory.DeclareInput for restriction mode requests A and P0
80 * PgPFactory.Build for prolongation mode calculates RowBasedOmega for prolongation operator
81 * PgPFactory.Build for restriction mode calculates RowBasedOmega for restriction operator
82 */
83 const ParameterList& pL = GetParameterList();
84 bool bReUseRowBasedOmegas = pL.get<bool>("ReUseRowBasedOmegas");
85 if (bReUseRowBasedOmegas == true && restrictionMode_ == true) {
86 coarseLevel.DeclareInput("RowBasedOmega", this, this); // RowBasedOmega is calculated by this PgPFactory and requested by this PgPFactory
87 }
88}
89
90template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
92 FactoryMonitor m(*this, "Prolongator smoothing (PG-AMG)", coarseLevel);
93
94 // Level Get
95 RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
96
97 // Get default tentative prolongator factory
98 // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
99 // -- Warning: Do not use directly initialPFact_. Use initialPFact instead everywhere!
100 RCP<const FactoryBase> initialPFact = GetFactory("P");
101 if (initialPFact == Teuchos::null) {
102 initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent");
103 }
104 RCP<Matrix> Ptent = coarseLevel.Get<RCP<Matrix> >("P", initialPFact.get());
105
107 if (restrictionMode_) {
108 SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
109 A = Utilities::Transpose(*A, true); // build transpose of A explicitely
110 }
111
113 bool doFillComplete = true;
114 bool optimizeStorage = true;
115 RCP<Matrix> DinvAP0 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A, false, *Ptent, false, GetOStream(Statistics2), doFillComplete, optimizeStorage);
116
117 doFillComplete = true;
118 optimizeStorage = false;
119 Teuchos::ArrayRCP<Scalar> diag = Utilities::GetMatrixDiagonal_arcp(*A);
120 Utilities::MyOldScaleMatrix(*DinvAP0, diag, true, doFillComplete, optimizeStorage); // scale matrix with reciprocal of diag
121
123
124 Teuchos::RCP<Vector> RowBasedOmega = Teuchos::null;
125
126 const ParameterList& pL = GetParameterList();
127 bool bReUseRowBasedOmegas = pL.get<bool>("ReUseRowBasedOmegas");
128 if (restrictionMode_ == false || bReUseRowBasedOmegas == false) {
129 // if in prolongation mode: calculate row based omegas
130 // if in restriction mode: calculate omegas only if row based omegas are not used from prolongation mode
131 ComputeRowBasedOmega(fineLevel, coarseLevel, A, Ptent, DinvAP0, RowBasedOmega);
132 } // if(bReUseRowBasedOmegas == false)
133 else {
134 // reuse row based omegas, calculated by this factory in the run before (with restrictionMode_ == false)
135 RowBasedOmega = coarseLevel.Get<Teuchos::RCP<Vector> >("RowBasedOmega", this);
136
137 // RowBasedOmega is now based on row map of A (not transposed)
138 // for restriction we use A^T instead of A
139 // -> recommunicate row based omega
140
141 // exporter: overlapping row map to nonoverlapping domain map (target map is unique)
142 // since A is already transposed we use the RangeMap of A
143 Teuchos::RCP<const Export> exporter =
144 ExportFactory::Build(RowBasedOmega->getMap(), A->getRangeMap());
145
146 Teuchos::RCP<Vector> noRowBasedOmega =
147 VectorFactory::Build(A->getRangeMap());
148
149 noRowBasedOmega->doExport(*RowBasedOmega, *exporter, Xpetra::INSERT);
150
151 // importer: nonoverlapping map to overlapping map
152
153 // importer: source -> target maps
154 Teuchos::RCP<const Import> importer =
155 ImportFactory::Build(A->getRangeMap(), A->getRowMap());
156
157 // doImport target->doImport(*source, importer, action)
158 RowBasedOmega->doImport(*noRowBasedOmega, *importer, Xpetra::INSERT);
159 }
160
161 Teuchos::ArrayRCP<Scalar> RowBasedOmega_local = RowBasedOmega->getDataNonConst(0);
162
164 RCP<Matrix> P_smoothed = Teuchos::null;
165 Utilities::MyOldScaleMatrix(*DinvAP0, RowBasedOmega_local, false, doFillComplete, optimizeStorage); // scale matrix with reciprocal of diag
166
167 Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TwoMatrixAdd(*Ptent, false, Teuchos::ScalarTraits<Scalar>::one(),
168 *DinvAP0, false, -Teuchos::ScalarTraits<Scalar>::one(),
169 P_smoothed, GetOStream(Statistics2));
170 P_smoothed->fillComplete(Ptent->getDomainMap(), Ptent->getRangeMap());
171
173
174 RCP<ParameterList> params = rcp(new ParameterList());
175 params->set("printLoadBalancingInfo", true);
176
177 // Level Set
178 if (!restrictionMode_) {
179 // prolongation factory is in prolongation mode
180 Set(coarseLevel, "P", P_smoothed);
181
182 // RfromPFactory used to indicate to TogglePFactory that a factory
183 // capable or producing R can be invoked later. TogglePFactory
184 // replaces dummy value with an index into it's array of prolongators
185 // pointing to the correct prolongator factory. This is later used by
186 // RfromP_Or_TransP to invoke the prolongatorfactory in RestrictionMode
187 int dummy = 7;
188 Set(coarseLevel, "RfromPfactory", dummy);
189
190 if (IsPrint(Statistics1))
191 GetOStream(Statistics1) << PerfUtils::PrintMatrixInfo(*P_smoothed, "P", params);
192
193 // NOTE: EXPERIMENTAL
194 if (Ptent->IsView("stridedMaps"))
195 P_smoothed->CreateView("stridedMaps", Ptent);
196
197 } else {
198 // prolongation factory is in restriction mode
199 RCP<Matrix> R = Utilities::Transpose(*P_smoothed, true);
200 Set(coarseLevel, "R", R);
201
202 if (IsPrint(Statistics1))
204
205 // NOTE: EXPERIMENTAL
206 if (Ptent->IsView("stridedMaps"))
207 R->CreateView("stridedMaps", Ptent, true);
208 }
209}
210
211template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
212void PgPFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::ComputeRowBasedOmega(Level& /* fineLevel */, Level& coarseLevel, const RCP<Matrix>& A, const RCP<Matrix>& P0, const RCP<Matrix>& DinvAP0, RCP<Vector>& RowBasedOmega) const {
213 FactoryMonitor m(*this, "PgPFactory::ComputeRowBasedOmega", coarseLevel);
214
215 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
216 Scalar sZero = Teuchos::ScalarTraits<Scalar>::zero();
217 Magnitude mZero = Teuchos::ScalarTraits<Scalar>::magnitude(sZero);
218
219 Teuchos::RCP<Vector> Numerator = Teuchos::null;
220 Teuchos::RCP<Vector> Denominator = Teuchos::null;
221
222 const ParameterList& pL = GetParameterList();
223 MinimizationNorm min_norm = pL.get<MinimizationNorm>("Minimization norm");
224
225 switch (min_norm) {
226 case ANORM: {
227 // MUEMAT mode (=paper)
228 // Minimize with respect to the (A)' A norm.
229 // Need to be smart here to avoid the construction of A' A
230 //
231 // diag( P0' (A' A) D^{-1} A P0)
232 // omega = ------------------------------------------
233 // diag( P0' A' D^{-1}' ( A' A) D^{-1} A P0)
234 //
235 // expensive, since we have to recalculate AP0 due to the lack of an explicit scaling routine for DinvAP0
236
237 // calculate A * P0
238 bool doFillComplete = true;
239 bool optimizeStorage = false;
240 RCP<Matrix> AP0 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A, false, *P0, false, GetOStream(Statistics2), doFillComplete, optimizeStorage);
241
242 // compute A * D^{-1} * A * P0
243 RCP<Matrix> ADinvAP0 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A, false, *DinvAP0, false, GetOStream(Statistics2), doFillComplete, optimizeStorage);
244
245 Numerator = VectorFactory::Build(ADinvAP0->getColMap(), true);
246 Denominator = VectorFactory::Build(ADinvAP0->getColMap(), true);
247 MultiplyAll(AP0, ADinvAP0, Numerator);
248 MultiplySelfAll(ADinvAP0, Denominator);
249 } break;
250 case L2NORM: {
251 // ML mode 1 (cheapest)
252 // Minimize with respect to L2 norm
253 // diag( P0' D^{-1} A P0)
254 // omega = -----------------------------
255 // diag( P0' A' D^{-1}' D^{-1} A P0)
256 //
257 Numerator = VectorFactory::Build(DinvAP0->getColMap(), true);
258 Denominator = VectorFactory::Build(DinvAP0->getColMap(), true);
259 MultiplyAll(P0, DinvAP0, Numerator);
260 MultiplySelfAll(DinvAP0, Denominator);
261 } break;
262 case DINVANORM: {
263 // ML mode 2
264 // Minimize with respect to the (D^{-1} A)' D^{-1} A norm.
265 // Need to be smart here to avoid the construction of A' A
266 //
267 // diag( P0' ( A' D^{-1}' D^{-1} A) D^{-1} A P0)
268 // omega = --------------------------------------------------------
269 // diag( P0' A' D^{-1}' ( A' D^{-1}' D^{-1} A) D^{-1} A P0)
270 //
271
272 // compute D^{-1} * A * D^{-1} * A * P0
273 bool doFillComplete = true;
274 bool optimizeStorage = true;
275 Teuchos::ArrayRCP<Scalar> diagA = Utilities::GetMatrixDiagonal_arcp(*A);
276 RCP<Matrix> DinvADinvAP0 = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*A, false, *DinvAP0, false, GetOStream(Statistics2), doFillComplete, optimizeStorage);
277 Utilities::MyOldScaleMatrix(*DinvADinvAP0, diagA, true, doFillComplete, optimizeStorage); // scale matrix with reciprocal of diag
278 diagA = Teuchos::ArrayRCP<Scalar>();
279
280 Numerator = VectorFactory::Build(DinvADinvAP0->getColMap(), true);
281 Denominator = VectorFactory::Build(DinvADinvAP0->getColMap(), true);
282 MultiplyAll(DinvAP0, DinvADinvAP0, Numerator);
283 MultiplySelfAll(DinvADinvAP0, Denominator);
284 } break;
285 default:
286 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::PgPFactory::Build: minimization mode not supported. error");
287 }
288
290 Teuchos::RCP<Vector> ColBasedOmega =
291 VectorFactory::Build(Numerator->getMap() /*DinvAP0->getColMap()*/, true);
292
293 ColBasedOmega->putScalar(-666 );
294
295 Teuchos::ArrayRCP<const Scalar> Numerator_local = Numerator->getData(0);
296 Teuchos::ArrayRCP<const Scalar> Denominator_local = Denominator->getData(0);
297 Teuchos::ArrayRCP<Scalar> ColBasedOmega_local = ColBasedOmega->getDataNonConst(0);
298 GlobalOrdinal zero_local = 0; // count negative colbased omegas
299 GlobalOrdinal nan_local = 0; // count NaNs -> set them to zero
300 Magnitude min_local = 1000000.0; // Teuchos::ScalarTraits<Scalar>::one() * (Scalar) 1000000;
301 Magnitude max_local = 0.0;
302 for (LocalOrdinal i = 0; i < Teuchos::as<LocalOrdinal>(Numerator->getLocalLength()); i++) {
303 if (Teuchos::ScalarTraits<Scalar>::magnitude(Denominator_local[i]) == mZero) {
304 ColBasedOmega_local[i] = 0.0; // fallback: nonsmoothed basis function since denominator == 0.0
305 nan_local++;
306 } else {
307 ColBasedOmega_local[i] = Numerator_local[i] / Denominator_local[i]; // default case
308 }
309
310 if (Teuchos::ScalarTraits<Scalar>::magnitude(ColBasedOmega_local[i]) < mZero) { // negative omegas are not valid. set them to zero
311 ColBasedOmega_local[i] = Teuchos::ScalarTraits<Scalar>::zero();
312 zero_local++; // count zero omegas
313 }
314
315 // handle case that Nominator == Denominator -> Dirichlet bcs in A?
316 // fallback if ColBasedOmega == 1 -> very strong smoothing may lead to zero rows in P
317 // TAW: this is somewhat nonstandard and a rough fallback strategy to avoid problems
318 // also avoid "overshooting" with omega > 0.8
319 if (Teuchos::ScalarTraits<Scalar>::magnitude(ColBasedOmega_local[i]) >= 0.8) {
320 ColBasedOmega_local[i] = 0.0;
321 }
322
323 if (Teuchos::ScalarTraits<Scalar>::magnitude(ColBasedOmega_local[i]) < min_local) {
324 min_local = Teuchos::ScalarTraits<Scalar>::magnitude(ColBasedOmega_local[i]);
325 }
326 if (Teuchos::ScalarTraits<Scalar>::magnitude(ColBasedOmega_local[i]) > max_local) {
327 max_local = Teuchos::ScalarTraits<Scalar>::magnitude(ColBasedOmega_local[i]);
328 }
329 }
330
331 { // be verbose
332 GlobalOrdinal zero_all;
333 GlobalOrdinal nan_all;
334 Magnitude min_all;
335 Magnitude max_all;
336 MueLu_sumAll(A->getRowMap()->getComm(), zero_local, zero_all);
337 MueLu_sumAll(A->getRowMap()->getComm(), nan_local, nan_all);
338 MueLu_minAll(A->getRowMap()->getComm(), min_local, min_all);
339 MueLu_maxAll(A->getRowMap()->getComm(), max_local, max_all);
340
341 GetOStream(MueLu::Statistics1, 0) << "PgPFactory: smoothed aggregation (scheme: ";
342 switch (min_norm) {
343 case ANORM: GetOStream(Statistics1) << "Anorm)" << std::endl; break;
344 case L2NORM: GetOStream(Statistics1) << "L2norm)" << std::endl; break;
345 case DINVANORM: GetOStream(Statistics1) << "DinvAnorm)" << std::endl; break;
346 default: GetOStream(Statistics1) << "unknown)" << std::endl; break;
347 }
348 GetOStream(Statistics1) << "Damping parameter: min = " << min_all << ", max = " << max_all << std::endl;
349 GetOStream(Statistics) << "# negative omegas: " << zero_all << " out of " << ColBasedOmega->getGlobalLength() << " column-based omegas" << std::endl;
350 GetOStream(Statistics) << "# NaNs: " << nan_all << " out of " << ColBasedOmega->getGlobalLength() << " column-based omegas" << std::endl;
351 }
352
353 if (coarseLevel.IsRequested("ColBasedOmega", this)) {
354 coarseLevel.Set("ColBasedOmega", ColBasedOmega, this);
355 }
356
358 // transform column based omegas to row based omegas
359 RowBasedOmega =
360 VectorFactory::Build(DinvAP0->getRowMap(), true);
361
362 RowBasedOmega->putScalar(-666); // TODO bad programming style
363
364 bool bAtLeastOneDefined = false;
365 Teuchos::ArrayRCP<Scalar> RowBasedOmega_local = RowBasedOmega->getDataNonConst(0);
366 for (LocalOrdinal row = 0; row < Teuchos::as<LocalOrdinal>(A->getLocalNumRows()); row++) {
367 Teuchos::ArrayView<const LocalOrdinal> lindices;
368 Teuchos::ArrayView<const Scalar> lvals;
369 DinvAP0->getLocalRowView(row, lindices, lvals);
370 bAtLeastOneDefined = false;
371 for (size_t j = 0; j < Teuchos::as<size_t>(lindices.size()); j++) {
372 Scalar omega = ColBasedOmega_local[lindices[j]];
373 if (Teuchos::ScalarTraits<Scalar>::magnitude(omega) != -666) { // TODO bad programming style
374 bAtLeastOneDefined = true;
375 if (Teuchos::ScalarTraits<Scalar>::magnitude(RowBasedOmega_local[row]) == -666)
376 RowBasedOmega_local[row] = omega;
377 else if (Teuchos::ScalarTraits<Scalar>::magnitude(omega) < Teuchos::ScalarTraits<Scalar>::magnitude(RowBasedOmega_local[row]))
378 RowBasedOmega_local[row] = omega;
379 }
380 }
381 if (bAtLeastOneDefined == true) {
382 if (Teuchos::ScalarTraits<Scalar>::magnitude(RowBasedOmega_local[row]) < mZero)
383 RowBasedOmega_local[row] = sZero;
384 }
385 }
386
387 if (coarseLevel.IsRequested("RowBasedOmega", this)) {
388 Set(coarseLevel, "RowBasedOmega", RowBasedOmega);
389 }
390}
391
392template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
393void PgPFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MultiplySelfAll(const RCP<Matrix>& Op, Teuchos::RCP<Vector>& InnerProdVec) const {
394 // note: InnerProdVec is based on column map of Op
395 TEUCHOS_TEST_FOR_EXCEPTION(!InnerProdVec->getMap()->isSameAs(*Op->getColMap()), Exceptions::RuntimeError, "MueLu::PgPFactory::MultiplySelfAll: map of InnerProdVec must be same as column map of operator. error");
396
397 Teuchos::ArrayRCP<Scalar> InnerProd_local = InnerProdVec->getDataNonConst(0);
398
399 Teuchos::ArrayView<const LocalOrdinal> lindices;
400 Teuchos::ArrayView<const Scalar> lvals;
401
402 for (size_t n = 0; n < Op->getLocalNumRows(); n++) {
403 Op->getLocalRowView(n, lindices, lvals);
404 for (size_t i = 0; i < Teuchos::as<size_t>(lindices.size()); i++) {
405 InnerProd_local[lindices[i]] += lvals[i] * lvals[i];
406 }
407 }
408 InnerProd_local = Teuchos::ArrayRCP<Scalar>();
409
410 // exporter: overlapping map to nonoverlapping map (target map is unique)
411 Teuchos::RCP<const Export> exporter =
412 ExportFactory::Build(Op->getColMap(), Op->getDomainMap());
413
414 Teuchos::RCP<Vector> nonoverlap =
415 VectorFactory::Build(Op->getDomainMap());
416
417 nonoverlap->doExport(*InnerProdVec, *exporter, Xpetra::ADD);
418
419 // importer: nonoverlapping map to overlapping map
420
421 // importer: source -> target maps
422 Teuchos::RCP<const Import> importer =
423 ImportFactory::Build(Op->getDomainMap(), Op->getColMap());
424
425 // doImport target->doImport(*source, importer, action)
426 InnerProdVec->doImport(*nonoverlap, *importer, Xpetra::INSERT);
427}
428
429template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
430void PgPFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::MultiplyAll(const RCP<Matrix>& left, const RCP<Matrix>& right, Teuchos::RCP<Vector>& InnerProdVec) const {
431 TEUCHOS_TEST_FOR_EXCEPTION(!left->getDomainMap()->isSameAs(*right->getDomainMap()), Exceptions::RuntimeError, "MueLu::PgPFactory::MultiplyAll: domain maps of left and right do not match. Error.");
432 TEUCHOS_TEST_FOR_EXCEPTION(!left->getRowMap()->isSameAs(*right->getRowMap()), Exceptions::RuntimeError, "MueLu::PgPFactory::MultiplyAll: row maps of left and right do not match. Error.");
433#if 1 // 1=new "fast code, 0=old "slow", but safe code
434#if 0 // not necessary - remove me
435 if(InnerProdVec->getMap()->isSameAs(*left->getColMap())) {
436 // initialize NewRightLocal vector and assign all entries to
437 // left->getColMap()->getLocalNumElements() + 1
438 std::vector<LocalOrdinal> NewRightLocal(right->getColMap()->getLocalNumElements(), Teuchos::as<LocalOrdinal>(left->getColMap()->getLocalNumElements()+1));
439
440 LocalOrdinal i = 0;
441 for (size_t j=0; j < right->getColMap()->getLocalNumElements(); j++) {
442 while ( (i < Teuchos::as<LocalOrdinal>(left->getColMap()->getLocalNumElements())) &&
443 (left->getColMap()->getGlobalElement(i) < right->getColMap()->getGlobalElement(j)) ) i++;
444 if (left->getColMap()->getGlobalElement(i) == right->getColMap()->getGlobalElement(j)) {
445 NewRightLocal[j] = i;
446 }
447 }
448
449 Teuchos::ArrayRCP< Scalar > InnerProd_local = InnerProdVec->getDataNonConst(0);
450 std::vector<Scalar> temp_array(left->getColMap()->getLocalNumElements()+1, 0.0);
451
452 for(size_t n=0; n<right->getLocalNumRows(); n++) {
453 Teuchos::ArrayView<const LocalOrdinal> lindices_left;
454 Teuchos::ArrayView<const Scalar> lvals_left;
455 Teuchos::ArrayView<const LocalOrdinal> lindices_right;
456 Teuchos::ArrayView<const Scalar> lvals_right;
457
458 left->getLocalRowView (n, lindices_left, lvals_left);
459 right->getLocalRowView(n, lindices_right, lvals_right);
460
461 for(size_t j=0; j<Teuchos::as<size_t>(lindices_right.size()); j++) {
462 temp_array[NewRightLocal[lindices_right[j] ] ] = lvals_right[j];
463 }
464 for (size_t j=0; j < Teuchos::as<size_t>(lindices_left.size()); j++) {
465 InnerProd_local[lindices_left[j]] += temp_array[lindices_left[j] ]*lvals_left[j];
466 }
467 for (size_t j=0; j < Teuchos::as<size_t>(lindices_right.size()); j++) {
468 temp_array[NewRightLocal[lindices_right[j] ] ] = 0.0;
469 }
470 }
471 // exporter: overlapping map to nonoverlapping map (target map is unique)
472 Teuchos::RCP<const Export> exporter =
473 ExportFactory::Build(left->getColMap(), left->getDomainMap()); // TODO: change left to right?
474
475 Teuchos::RCP<Vector > nonoverlap =
476 VectorFactory::Build(left->getDomainMap()); // TODO: change left to right?
477
478 nonoverlap->doExport(*InnerProdVec, *exporter, Xpetra::ADD);
479
480 // importer: nonoverlapping map to overlapping map
481
482 // importer: source -> target maps
483 Teuchos::RCP<const Import > importer =
484 ImportFactory::Build(left->getDomainMap(), left->getColMap()); // TODO: change left to right?
485
486 // doImport target->doImport(*source, importer, action)
487 InnerProdVec->doImport(*nonoverlap, *importer, Xpetra::INSERT);
488
489
490 } else
491#endif // end remove me
492 if (InnerProdVec->getMap()->isSameAs(*right->getColMap())) {
493 size_t szNewLeftLocal = TEUCHOS_MAX(left->getColMap()->getLocalNumElements(), right->getColMap()->getLocalNumElements());
494 Teuchos::RCP<std::vector<LocalOrdinal> > NewLeftLocal = Teuchos::rcp(new std::vector<LocalOrdinal>(szNewLeftLocal, Teuchos::as<LocalOrdinal>(right->getColMap()->getMaxLocalIndex() + 1)));
495
496 LocalOrdinal j = 0;
497 for (size_t i = 0; i < left->getColMap()->getLocalNumElements(); i++) {
498 while ((j < Teuchos::as<LocalOrdinal>(right->getColMap()->getLocalNumElements())) &&
499 (right->getColMap()->getGlobalElement(j) < left->getColMap()->getGlobalElement(i))) j++;
500 if (right->getColMap()->getGlobalElement(j) == left->getColMap()->getGlobalElement(i)) {
501 (*NewLeftLocal)[i] = j;
502 }
503 }
504
505 /*for (size_t i=0; i < right->getColMap()->getLocalNumElements(); i++) {
506 std::cout << "left col map: " << (*NewLeftLocal)[i] << " GID: " << left->getColMap()->getGlobalElement((*NewLeftLocal)[i]) << " GID: " << right->getColMap()->getGlobalElement(i) << " right col map: " << i << std::endl;
507 }*/
508
509 Teuchos::ArrayRCP<Scalar> InnerProd_local = InnerProdVec->getDataNonConst(0);
510 Teuchos::RCP<std::vector<Scalar> > temp_array = Teuchos::rcp(new std::vector<Scalar>(right->getColMap()->getMaxLocalIndex() + 2, 0.0));
511
512 for (size_t n = 0; n < left->getLocalNumRows(); n++) {
513 Teuchos::ArrayView<const LocalOrdinal> lindices_left;
514 Teuchos::ArrayView<const Scalar> lvals_left;
515 Teuchos::ArrayView<const LocalOrdinal> lindices_right;
516 Teuchos::ArrayView<const Scalar> lvals_right;
517
518 left->getLocalRowView(n, lindices_left, lvals_left);
519 right->getLocalRowView(n, lindices_right, lvals_right);
520
521 for (size_t i = 0; i < Teuchos::as<size_t>(lindices_left.size()); i++) {
522 (*temp_array)[(*NewLeftLocal)[lindices_left[i]]] = lvals_left[i];
523 }
524 for (size_t i = 0; i < Teuchos::as<size_t>(lindices_right.size()); i++) {
525 InnerProd_local[lindices_right[i]] += (*temp_array)[lindices_right[i]] * lvals_right[i];
526 }
527 for (size_t i = 0; i < Teuchos::as<size_t>(lindices_left.size()); i++) {
528 (*temp_array)[(*NewLeftLocal)[lindices_left[i]]] = 0.0;
529 }
530 }
531 InnerProd_local = Teuchos::ArrayRCP<Scalar>();
532 // exporter: overlapping map to nonoverlapping map (target map is unique)
533 Teuchos::RCP<const Export> exporter =
534 ExportFactory::Build(right->getColMap(), right->getDomainMap()); // TODO: change left to right?
535
536 Teuchos::RCP<Vector> nonoverlap =
537 VectorFactory::Build(right->getDomainMap()); // TODO: change left to right?
538
539 nonoverlap->doExport(*InnerProdVec, *exporter, Xpetra::ADD);
540
541 // importer: nonoverlapping map to overlapping map
542
543 // importer: source -> target maps
544 Teuchos::RCP<const Import> importer =
545 ImportFactory::Build(right->getDomainMap(), right->getColMap()); // TODO: change left to right?
546 // doImport target->doImport(*source, importer, action)
547 InnerProdVec->doImport(*nonoverlap, *importer, Xpetra::INSERT);
548 } else {
549 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::PgPFactory::MultiplyAll: map of InnerProdVec must be same as column map of left operator? Error.");
550 }
551
552#else // old "safe" code
553 if (InnerProdVec->getMap()->isSameAs(*left->getColMap())) {
554 Teuchos::ArrayRCP<Scalar> InnerProd_local = InnerProdVec->getDataNonConst(0);
555
556 // declare variables
557 Teuchos::ArrayView<const LocalOrdinal> lindices_left;
558 Teuchos::ArrayView<const Scalar> lvals_left;
559 Teuchos::ArrayView<const LocalOrdinal> lindices_right;
560 Teuchos::ArrayView<const Scalar> lvals_right;
561
562 for (size_t n = 0; n < left->getLocalNumRows(); n++) {
563 left->getLocalRowView(n, lindices_left, lvals_left);
564 right->getLocalRowView(n, lindices_right, lvals_right);
565
566 for (size_t i = 0; i < Teuchos::as<size_t>(lindices_left.size()); i++) {
567 GlobalOrdinal left_gid = left->getColMap()->getGlobalElement(lindices_left[i]);
568 for (size_t j = 0; j < Teuchos::as<size_t>(lindices_right.size()); j++) {
569 GlobalOrdinal right_gid = right->getColMap()->getGlobalElement(lindices_right[j]);
570 if (left_gid == right_gid) {
571 InnerProd_local[lindices_left[i]] += lvals_left[i] * lvals_right[j];
572 break; // skip remaining gids of right operator
573 }
574 }
575 }
576 }
577
578 // exporter: overlapping map to nonoverlapping map (target map is unique)
579 Teuchos::RCP<const Export> exporter =
580 ExportFactory::Build(left->getColMap(), left->getDomainMap()); // TODO: change left to right?
581
582 Teuchos::RCP<Vector> nonoverlap =
583 VectorFactory::Build(left->getDomainMap()); // TODO: change left to right?
584
585 nonoverlap->doExport(*InnerProdVec, *exporter, Xpetra::ADD);
586
587 // importer: nonoverlapping map to overlapping map
588
589 // importer: source -> target maps
590 Teuchos::RCP<const Import> importer =
591 ImportFactory::Build(left->getDomainMap(), left->getColMap()); // TODO: change left to right?
592
593 // doImport target->doImport(*source, importer, action)
594 InnerProdVec->doImport(*nonoverlap, *importer, Xpetra::INSERT);
595 } else if (InnerProdVec->getMap()->isSameAs(*right->getColMap())) {
596 Teuchos::ArrayRCP<Scalar> InnerProd_local = InnerProdVec->getDataNonConst(0);
597
598 Teuchos::ArrayView<const LocalOrdinal> lindices_left;
599 Teuchos::ArrayView<const Scalar> lvals_left;
600 Teuchos::ArrayView<const LocalOrdinal> lindices_right;
601 Teuchos::ArrayView<const Scalar> lvals_right;
602
603 for (size_t n = 0; n < left->getLocalNumRows(); n++) {
604 left->getLocalRowView(n, lindices_left, lvals_left);
605 right->getLocalRowView(n, lindices_right, lvals_right);
606
607 for (size_t i = 0; i < Teuchos::as<size_t>(lindices_left.size()); i++) {
608 GlobalOrdinal left_gid = left->getColMap()->getGlobalElement(lindices_left[i]);
609 for (size_t j = 0; j < Teuchos::as<size_t>(lindices_right.size()); j++) {
610 GlobalOrdinal right_gid = right->getColMap()->getGlobalElement(lindices_right[j]);
611 if (left_gid == right_gid) {
612 InnerProd_local[lindices_right[j]] += lvals_left[i] * lvals_right[j];
613 break; // skip remaining gids of right operator
614 }
615 }
616 }
617 }
618
619 // exporter: overlapping map to nonoverlapping map (target map is unique)
620 Teuchos::RCP<const Export> exporter =
621 ExportFactory::Build(right->getColMap(), right->getDomainMap()); // TODO: change left to right?
622
623 Teuchos::RCP<Vector> nonoverlap =
624 VectorFactory::Build(right->getDomainMap()); // TODO: change left to right?
625
626 nonoverlap->doExport(*InnerProdVec, *exporter, Xpetra::ADD);
627
628 // importer: nonoverlapping map to overlapping map
629
630 // importer: source -> target maps
631 Teuchos::RCP<const Import> importer =
632 ImportFactory::Build(right->getDomainMap(), right->getColMap()); // TODO: change left to right?
633
634 // doImport target->doImport(*source, importer, action)
635 InnerProdVec->doImport(*nonoverlap, *importer, Xpetra::INSERT);
636 } else {
637 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "MueLu::PgPFactory::MultiplyAll: map of InnerProdVec must be same as column map of left or right operator? Error.");
638 }
639#endif
640}
641
642template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
643void PgPFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildP(Level& /* fineLevel */, Level& /* coarseLevel */) const {
644 std::cout << "TODO: remove me" << std::endl;
645}
646
647template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
649 SetParameter("ReUseRowBasedOmegas", ParameterEntry(bReuse));
650}
651
652} // namespace MueLu
653
654#endif /* MUELU_PGPFACTORY_DEF_HPP */
#define MueLu_maxAll(rcpComm, in, out)
#define MueLu_sumAll(rcpComm, in, out)
#define MueLu_minAll(rcpComm, in, out)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
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.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
const RCP< const FactoryBase > GetFactory(const std::string &varName) const
Default implementation of FactoryAcceptor::GetFactory().
Class that holds all level-specific information.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput().
const RCP< const FactoryManagerBase > GetFactoryManager()
returns the current factory manager
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
bool IsRequested(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need has been requested. Note: this tells nothing about whether the need's value exist...
virtual const Teuchos::ParameterList & GetParameterList() const
void SetParameter(const std::string &name, const ParameterEntry &entry)
Set a parameter directly as a ParameterEntry.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void ReUseDampingParameters(bool bReuse)
void MultiplyAll(const RCP< Matrix > &left, const RCP< Matrix > &right, Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &InnerProdVec) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void ComputeRowBasedOmega(Level &fineLevel, Level &coarseLevel, const RCP< Matrix > &A, const RCP< Matrix > &P0, const RCP< Matrix > &DinvAP0, RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &RowBasedOmega) const
void SetMinimizationMode(MinimizationNorm minnorm)
Set minimization mode (L2NORM for cheapest, ANORM more expensive, DINVANORM = default).
MinimizationNorm GetMinimizationMode()
return minimization mode
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void MultiplySelfAll(const RCP< Matrix > &Op, Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &InnerProdVec) const
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
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)
static void MyOldScaleMatrix(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Teuchos::ArrayRCP< const Scalar > &scalingVector, bool doInverse=true, bool doFillComplete=true, bool doOptimizeStorage=true)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Statistics1
Print more statistics.
@ Statistics
Print all statistics.