MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Hierarchy_def.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// MueLu: A package for multigrid based preconditioning
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact
39// Jonathan Hu (jhu@sandia.gov)
40// Andrey Prokopenko (aprokop@sandia.gov)
41// Ray Tuminaro (rstumin@sandia.gov)
42//
43// ***********************************************************************
44//
45// @HEADER
46
47#ifndef MUELU_HIERARCHY_DEF_HPP
48#define MUELU_HIERARCHY_DEF_HPP
49
50#include <time.h>
51
52#include <algorithm>
53#include <sstream>
54
55#include <Xpetra_Matrix.hpp>
56#include <Xpetra_MultiVectorFactory.hpp>
57#include <Xpetra_Operator.hpp>
58#include <Xpetra_IO.hpp>
59
61
62#include "MueLu_FactoryManager.hpp"
63#include "MueLu_HierarchyUtils.hpp"
64#include "MueLu_TopRAPFactory.hpp"
65#include "MueLu_TopSmootherFactory.hpp"
66#include "MueLu_Level.hpp"
67#include "MueLu_Monitor.hpp"
68#include "MueLu_PerfUtils.hpp"
69#include "MueLu_PFactory.hpp"
70#include "MueLu_SmootherFactory.hpp"
72
73#include "Teuchos_TimeMonitor.hpp"
74
75
76
77namespace MueLu {
78
79 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
89
90 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
92 : Hierarchy()
93 {
94 setObjectLabel(label);
95 Levels_[0]->setObjectLabel(label);
96 }
97
98 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
113
114 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115 Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Hierarchy(const RCP<Matrix>& A, const std::string& label)
116 : Hierarchy(A)
117 {
118 setObjectLabel(label);
119 Levels_[0]->setObjectLabel(label);
120 }
121
122 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
124 int levelID = LastLevelID() + 1; // ID of the inserted level
125
126 if (level->GetLevelID() != -1 && (level->GetLevelID() != levelID))
127 GetOStream(Warnings1) << "Hierarchy::AddLevel(): Level with ID=" << level->GetLevelID() <<
128 " have been added at the end of the hierarchy\n but its ID have been redefined" <<
129 " because last level ID of the hierarchy was " << LastLevelID() << "." << std::endl;
130
131 Levels_.push_back(level);
132 level->SetLevelID(levelID);
133 level->setlib(lib_);
134
135 level->SetPreviousLevel( (levelID == 0) ? Teuchos::null : Levels_[LastLevelID() - 1] );
136 level->setObjectLabel(this->getObjectLabel());
137 }
138
139 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
141 RCP<Level> newLevel = Levels_[LastLevelID()]->Build(); // new coarse level, using copy constructor
142 newLevel->setlib(lib_);
143 this->AddLevel(newLevel); // add to hierarchy
144 }
145
146 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148 TEUCHOS_TEST_FOR_EXCEPTION(levelID < 0 || levelID > LastLevelID(), Exceptions::RuntimeError,
149 "MueLu::Hierarchy::GetLevel(): invalid input parameter value: LevelID = " << levelID);
150 return Levels_[levelID];
151 }
152
153 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
157
158 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
160 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator> >("A");
161 RCP<const Teuchos::Comm<int> > comm = A->getDomainMap()->getComm();
162
163 int numLevels = GetNumLevels();
164 int numGlobalLevels;
165 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, numLevels, Teuchos::ptr(&numGlobalLevels));
166
167 return numGlobalLevels;
168 }
169
170 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
172 double totalNnz = 0, lev0Nnz = 1;
173 for (int i = 0; i < GetNumLevels(); ++i) {
174 TEUCHOS_TEST_FOR_EXCEPTION(!(Levels_[i]->IsAvailable("A")) , Exceptions::RuntimeError,
175 "Operator complexity cannot be calculated because A is unavailable on level " << i);
176 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator> >("A");
177 if (A.is_null())
178 break;
179
180 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
181 if (Am.is_null()) {
182 GetOStream(Warnings0) << "Some level operators are not matrices, operator complexity calculation aborted" << std::endl;
183 return 0.0;
184 }
185
186 totalNnz += as<double>(Am->getGlobalNumEntries());
187 if (i == 0)
188 lev0Nnz = totalNnz;
189 }
190 return totalNnz / lev0Nnz;
191 }
192
193 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
195 double node_sc = 0, global_sc=0;
196 double a0_nnz =0;
197 const size_t INVALID = Teuchos::OrdinalTraits<size_t>::invalid();
198 // Get cost of fine matvec
199 if (GetNumLevels() <= 0) return -1.0;
200 if (!Levels_[0]->IsAvailable("A")) return -1.0;
201
202 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator> >("A");
203 if (A.is_null()) return -1.0;
204 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
205 if(Am.is_null()) return -1.0;
206 a0_nnz = as<double>(Am->getGlobalNumEntries());
207
208 // Get smoother complexity at each level
209 for (int i = 0; i < GetNumLevels(); ++i) {
210 size_t level_sc=0;
211 if(!Levels_[i]->IsAvailable("PreSmoother")) continue;
212 RCP<SmootherBase> S = Levels_[i]->template Get<RCP<SmootherBase> >("PreSmoother");
213 if (S.is_null()) continue;
214 level_sc = S->getNodeSmootherComplexity();
215 if(level_sc == INVALID) {global_sc=-1.0;break;}
216
217 node_sc += as<double>(level_sc);
218 }
219
220 double min_sc=0.0;
221 RCP<const Teuchos::Comm<int> > comm =A->getDomainMap()->getComm();
222 Teuchos::reduceAll(*comm,Teuchos::REDUCE_SUM,node_sc,Teuchos::ptr(&global_sc));
223 Teuchos::reduceAll(*comm,Teuchos::REDUCE_MIN,node_sc,Teuchos::ptr(&min_sc));
224
225 if(min_sc < 0.0) return -1.0;
226 else return global_sc / a0_nnz;
227 }
228
229
230
231
232 // Coherence checks todo in Setup() (using an helper function):
233 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
235 TEUCHOS_TEST_FOR_EXCEPTION(level.lib() != lib_, Exceptions::RuntimeError,
236 "MueLu::Hierarchy::CheckLevel(): wrong underlying linear algebra library.");
237 TEUCHOS_TEST_FOR_EXCEPTION(level.GetLevelID() != levelID, Exceptions::RuntimeError,
238 "MueLu::Hierarchy::CheckLevel(): wrong level ID");
239 TEUCHOS_TEST_FOR_EXCEPTION(levelID != 0 && level.GetPreviousLevel() != Levels_[levelID-1], Exceptions::RuntimeError,
240 "MueLu::Hierarchy::Setup(): wrong level parent");
241 }
242
243 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
245 for (int i = 0; i < GetNumLevels(); ++i) {
246 RCP<Level> level = Levels_[i];
247 if (level->IsAvailable("A")) {
248 RCP<Operator> Aop = level->Get<RCP<Operator> >("A");
249 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Aop);
250 if (!A.is_null()) {
251 RCP<const Import> xpImporter = A->getCrsGraph()->getImporter();
252 if (!xpImporter.is_null())
253 xpImporter->setDistributorParameters(matvecParams);
254 RCP<const Export> xpExporter = A->getCrsGraph()->getExporter();
255 if (!xpExporter.is_null())
256 xpExporter->setDistributorParameters(matvecParams);
257 }
258 }
259 if (level->IsAvailable("P")) {
260 RCP<Matrix> P = level->Get<RCP<Matrix> >("P");
261 RCP<const Import> xpImporter = P->getCrsGraph()->getImporter();
262 if (!xpImporter.is_null())
263 xpImporter->setDistributorParameters(matvecParams);
264 RCP<const Export> xpExporter = P->getCrsGraph()->getExporter();
265 if (!xpExporter.is_null())
266 xpExporter->setDistributorParameters(matvecParams);
267 }
268 if (level->IsAvailable("R")) {
269 RCP<Matrix> R = level->Get<RCP<Matrix> >("R");
270 RCP<const Import> xpImporter = R->getCrsGraph()->getImporter();
271 if (!xpImporter.is_null())
272 xpImporter->setDistributorParameters(matvecParams);
273 RCP<const Export> xpExporter = R->getCrsGraph()->getExporter();
274 if (!xpExporter.is_null())
275 xpExporter->setDistributorParameters(matvecParams);
276 }
277 if (level->IsAvailable("Importer")) {
278 RCP<const Import> xpImporter = level->Get< RCP<const Import> >("Importer");
279 if (!xpImporter.is_null())
280 xpImporter->setDistributorParameters(matvecParams);
281 }
282 }
283 }
284
285 // The function uses three managers: fine, coarse and next coarse
286 // We construct the data for the coarse level, and do requests for the next coarse
287 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289 const RCP<const FactoryManagerBase> fineLevelManager,
290 const RCP<const FactoryManagerBase> coarseLevelManager,
291 const RCP<const FactoryManagerBase> nextLevelManager) {
292 // Use PrintMonitor/TimerMonitor instead of just a FactoryMonitor to print "Level 0" instead of Hierarchy(0)
293 // Print is done after the requests for next coarse level
294
295 TEUCHOS_TEST_FOR_EXCEPTION(LastLevelID() < coarseLevelID, Exceptions::RuntimeError,
296 "MueLu::Hierarchy:Setup(): level " << coarseLevelID << " (specified by coarseLevelID argument) "
297 "must be built before calling this function.");
298
299 Level& level = *Levels_[coarseLevelID];
300
301 std::string label = FormattingHelper::getColonLabel(level.getObjectLabel());
302 TimeMonitor m1(*this, label + this->ShortClassName() + ": " + "Setup (total)");
303 TimeMonitor m2(*this, label + this->ShortClassName() + ": " + "Setup" + " (total, level=" + Teuchos::toString(coarseLevelID) + ")");
304
305 // TODO: pass coarseLevelManager by reference
306 TEUCHOS_TEST_FOR_EXCEPTION(coarseLevelManager == Teuchos::null, Exceptions::RuntimeError,
307 "MueLu::Hierarchy::Setup(): argument coarseLevelManager cannot be null");
308
311
312 if (levelManagers_.size() < coarseLevelID+1)
313 levelManagers_.resize(coarseLevelID+1);
314 levelManagers_[coarseLevelID] = coarseLevelManager;
315
316 bool isFinestLevel = (fineLevelManager.is_null());
317 bool isLastLevel = (nextLevelManager.is_null());
318
319 int oldRank = -1;
320 if (isFinestLevel) {
321 RCP<Operator> A = level.Get< RCP<Operator> >("A");
322 RCP<const Map> domainMap = A->getDomainMap();
323 RCP<const Teuchos::Comm<int> > comm = domainMap->getComm();
324
325 // Initialize random seed for reproducibility
327
328 // Record the communicator on the level (used for timers sync)
329 level.SetComm(comm);
330 oldRank = SetProcRankVerbose(comm->getRank());
331
332 // Set the Hierarchy library to match that of the finest level matrix,
333 // even if it was already set
334 lib_ = domainMap->lib();
335 level.setlib(lib_);
336
337 } else {
338 // Permeate library to a coarser level
339 level.setlib(lib_);
340
341 Level& prevLevel = *Levels_[coarseLevelID-1];
342 oldRank = SetProcRankVerbose(prevLevel.GetComm()->getRank());
343 }
344
345 CheckLevel(level, coarseLevelID);
346
347 // Attach FactoryManager to the fine level
348 RCP<SetFactoryManager> SFMFine;
349 if (!isFinestLevel)
350 SFMFine = rcp(new SetFactoryManager(Levels_[coarseLevelID-1], fineLevelManager));
351
352 if (isFinestLevel && Levels_[coarseLevelID]->IsAvailable("Coordinates"))
353 ReplaceCoordinateMap(*Levels_[coarseLevelID]);
354
355 // Attach FactoryManager to the coarse level
356 SetFactoryManager SFMCoarse(Levels_[coarseLevelID], coarseLevelManager);
357
358 if (isDumpingEnabled_ && (dumpLevel_ == 0 || dumpLevel_ == -1) && coarseLevelID == 1)
360
361 RCP<TopSmootherFactory> coarseFact;
362 RCP<TopSmootherFactory> smootherFact = rcp(new TopSmootherFactory(coarseLevelManager, "Smoother"));
363
364 int nextLevelID = coarseLevelID + 1;
365
366 RCP<SetFactoryManager> SFMNext;
367 if (isLastLevel == false) {
368 // We are not at the coarsest level, so there is going to be another level ("next coarse") after this one ("coarse")
369 if (nextLevelID > LastLevelID())
370 AddNewLevel();
371 CheckLevel(*Levels_[nextLevelID], nextLevelID);
372
373 // Attach FactoryManager to the next level (level after coarse)
374 SFMNext = rcp(new SetFactoryManager(Levels_[nextLevelID], nextLevelManager));
375 Levels_[nextLevelID]->Request(TopRAPFactory(coarseLevelManager, nextLevelManager));
376
377 // Do smoother requests here. We don't know whether this is going to be
378 // the coarsest level or not, but we need to DeclareInput before we call
379 // coarseRAPFactory.Build(), otherwise some stuff may be erased after
380 // level releases
381 level.Request(*smootherFact);
382
383 } else {
384 // Similar to smoother above, do the coarse solver request here. We don't
385 // know whether this is going to be the coarsest level or not, but we
386 // need to DeclareInput before we call coarseRAPFactory.Build(),
387 // otherwise some stuff may be erased after level releases. This is
388 // actually evident on ProjectorSmoother. It requires both "A" and
389 // "Nullspace". However, "Nullspace" is erased after all releases, so if
390 // we call the coarse factory request after RAP build we would not have
391 // any data, and cannot get it as we don't have previous managers. The
392 // typical trace looks like this:
393 //
394 // MueLu::Level(0)::GetFactory(Aggregates, 0): No FactoryManager
395 // during request for data " Aggregates" on level 0 by factory TentativePFactory
396 // during request for data " P" on level 1 by factory EminPFactory
397 // during request for data " P" on level 1 by factory TransPFactory
398 // during request for data " R" on level 1 by factory RAPFactory
399 // during request for data " A" on level 1 by factory TentativePFactory
400 // during request for data " Nullspace" on level 2 by factory NullspaceFactory
401 // during request for data " Nullspace" on level 2 by factory NullspacePresmoothFactory
402 // during request for data " Nullspace" on level 2 by factory ProjectorSmoother
403 // during request for data " PreSmoother" on level 2 by factory NoFactory
404 if (coarseFact.is_null())
405 coarseFact = rcp(new TopSmootherFactory(coarseLevelManager, "CoarseSolver"));
406 level.Request(*coarseFact);
407 }
408
409 GetOStream(Runtime0) << std::endl;
410 PrintMonitor m0(*this, "Level " + Teuchos::toString(coarseLevelID), static_cast<MsgType>(Runtime0 | Test));
411
412 // Build coarse level hierarchy
413 RCP<Operator> Ac = Teuchos::null;
414 TopRAPFactory coarseRAPFactory(fineLevelManager, coarseLevelManager);
415
416 if (level.IsAvailable("A")) {
417 Ac = level.Get<RCP<Operator> >("A");
418 } else if (!isFinestLevel) {
419 // We only build here, the release is done later
420 coarseRAPFactory.Build(*level.GetPreviousLevel(), level);
421 }
422
423 bool setLastLevelviaMaxCoarseSize = false;
424 if (level.IsAvailable("A"))
425 Ac = level.Get<RCP<Operator> >("A");
426 RCP<Matrix> Acm = rcp_dynamic_cast<Matrix>(Ac);
427
428 // Record the communicator on the level
429 if (!Ac.is_null())
430 level.SetComm(Ac->getDomainMap()->getComm());
431
432 // Test if we reach the end of the hierarchy
433 bool isOrigLastLevel = isLastLevel;
434 if (isLastLevel) {
435 // Last level as we have achieved the max limit
436 isLastLevel = true;
437
438 } else if (Ac.is_null()) {
439 // Last level for this processor, as it does not belong to the next
440 // subcommunicator. Other processors may continue working on the
441 // hierarchy
442 isLastLevel = true;
443
444 } else {
445 if (!Acm.is_null() && Acm->getGlobalNumRows() <= maxCoarseSize_) {
446 // Last level as the size of the coarse matrix became too small
447 GetOStream(Runtime0) << "Max coarse size (<= " << maxCoarseSize_ << ") achieved" << std::endl;
448 isLastLevel = true;
449 if (Acm->getGlobalNumRows() != 0) setLastLevelviaMaxCoarseSize = true;
450 }
451 }
452
453 if (!Ac.is_null() && !isFinestLevel) {
454 RCP<Operator> A = Levels_[coarseLevelID-1]->template Get< RCP<Operator> >("A");
455 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
456
457 const double maxCoarse2FineRatio = 0.8;
458 if (!Acm.is_null() && !Am.is_null() && Acm->getGlobalNumRows() > maxCoarse2FineRatio * Am->getGlobalNumRows()) {
459 // We could abort here, but for now we simply notify user.
460 // Couple of additional points:
461 // - if repartitioning is delayed until level K, but the aggregation
462 // procedure stagnates between levels K-1 and K. In this case,
463 // repartitioning could enable faster coarsening once again, but the
464 // hierarchy construction will abort due to the stagnation check.
465 // - if the matrix is small enough, we could move it to one processor.
466 GetOStream(Warnings0) << "Aggregation stagnated. Please check your matrix and/or adjust your configuration file."
467 << "Possible fixes:\n"
468 << " - reduce the maximum number of levels\n"
469 << " - enable repartitioning\n"
470 << " - increase the minimum coarse size." << std::endl;
471
472 }
473 }
474
475 if (isLastLevel) {
476 if (!isOrigLastLevel) {
477 // We did not expect to finish this early so we did request a smoother.
478 // We need a coarse solver instead. Do the magic.
479 level.Release(*smootherFact);
480 if (coarseFact.is_null())
481 coarseFact = rcp(new TopSmootherFactory(coarseLevelManager, "CoarseSolver"));
482 level.Request(*coarseFact);
483 }
484
485 // Do the actual build, if we have any data.
486 // NOTE: this is not a great check, we may want to call Build() regardless.
487 if (!Ac.is_null())
488 coarseFact->Build(level);
489
490 // Once the dirty deed is done, release stuff. The smoother has already
491 // been released.
492 level.Release(*coarseFact);
493
494 } else {
495 // isLastLevel = false => isOrigLastLevel = false, meaning that we have
496 // requested the smoother. Now we need to build it and to release it.
497 // We don't need to worry about the coarse solver, as we didn't request it.
498 if (!Ac.is_null())
499 smootherFact->Build(level);
500
501 level.Release(*smootherFact);
502 }
503
504 if (isLastLevel == true) {
505 int actualNumLevels = nextLevelID;
506 if (isOrigLastLevel == false) {
507 // Earlier in the function, we constructed the next coarse level, and requested data for the that level,
508 // assuming that we are not at the coarsest level. Now, we changed our mind, so we have to release those.
509 Levels_[nextLevelID]->Release(TopRAPFactory(coarseLevelManager, nextLevelManager));
510
511 // We truncate/resize the hierarchy and possibly remove the last created level if there is
512 // something wrong with it as indicated by its P not being valid. This might happen
513 // if the global number of aggregates turns out to be zero
514
515
516 if (!setLastLevelviaMaxCoarseSize) {
517 if (Levels_[nextLevelID-1]->IsAvailable("P")) {
518 if (Levels_[nextLevelID-1]->template Get<RCP<Matrix> >("P") == Teuchos::null) actualNumLevels = nextLevelID-1;
519 }
520 else actualNumLevels = nextLevelID-1;
521 }
522 }
523 if (actualNumLevels == nextLevelID-1) {
524 // Didn't expect to finish early so we requested smoother but need coarse solver instead.
525 Levels_[nextLevelID-2]->Release(*smootherFact);
526
527 if (Levels_[nextLevelID-2]->IsAvailable("PreSmoother") ) Levels_[nextLevelID-2]->RemoveKeepFlag("PreSmoother" ,NoFactory::get());
528 if (Levels_[nextLevelID-2]->IsAvailable("PostSmoother")) Levels_[nextLevelID-2]->RemoveKeepFlag("PostSmoother",NoFactory::get());
529 if (coarseFact.is_null())
530 coarseFact = rcp(new TopSmootherFactory(coarseLevelManager, "CoarseSolver"));
531 Levels_[nextLevelID-2]->Request(*coarseFact);
532 if ( !(Levels_[nextLevelID-2]->template Get<RCP<Matrix> >("A").is_null() ))
533 coarseFact->Build( *(Levels_[nextLevelID-2]));
534 Levels_[nextLevelID-2]->Release(*coarseFact);
535 }
536 Levels_.resize(actualNumLevels);
537 }
538
539 // I think this is the proper place for graph so that it shows every dependence
540 if (isDumpingEnabled_ && ( (dumpLevel_ > 0 && coarseLevelID == dumpLevel_) || dumpLevel_ == -1 ) )
541 DumpCurrentGraph(coarseLevelID);
542
543 if (!isFinestLevel) {
544 // Release the hierarchy data
545 // We release so late to help blocked solvers, as the smoothers for them need A blocks
546 // which we construct in RAPFactory
547 level.Release(coarseRAPFactory);
548 }
549
550 if (oldRank != -1)
551 SetProcRankVerbose(oldRank);
552
553 return isLastLevel;
554 }
555
556 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
558 int numLevels = Levels_.size();
559 TEUCHOS_TEST_FOR_EXCEPTION(levelManagers_.size() != numLevels, Exceptions::RuntimeError,
560 "Hierarchy::SetupRe: " << Levels_.size() << " levels, but " << levelManagers_.size() << " level factory managers");
561
562 const int startLevel = 0;
563 Clear(startLevel);
564
565#ifdef HAVE_MUELU_DEBUG
566 // Reset factories' data used for debugging
567 for (int i = 0; i < numLevels; i++)
568 levelManagers_[i]->ResetDebugData();
569
570#endif
571
572 int levelID;
573 for (levelID = startLevel; levelID < numLevels;) {
574 bool r = Setup(levelID,
575 (levelID != 0 ? levelManagers_[levelID-1] : Teuchos::null),
576 levelManagers_[levelID],
577 (levelID+1 != numLevels ? levelManagers_[levelID+1] : Teuchos::null));
578 levelID++;
579 if (r) break;
580 }
581 // We may construct fewer levels for some reason, make sure we continue
582 // doing that in the future
583 Levels_ .resize(levelID);
584 levelManagers_.resize(levelID);
585
586 int sizeOfVecs = sizeOfAllocatedLevelMultiVectors_;
587
588 AllocateLevelMultiVectors(sizeOfVecs, true);
589
590 // since the # of levels, etc. may have changed, force re-determination of description during next call to description()
592
594 }
595
596 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
597 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Setup(const FactoryManagerBase& manager, int startLevel, int numDesiredLevels) {
598 // Use MueLu::BaseClass::description() to avoid printing "{numLevels = 1}" (numLevels is increasing...)
599 PrintMonitor m0(*this, "Setup (" + this->MueLu::BaseClass::description() + ")", Runtime0);
600
601 Clear(startLevel);
602
603 // Check Levels_[startLevel] exists.
604 TEUCHOS_TEST_FOR_EXCEPTION(Levels_.size() <= startLevel, Exceptions::RuntimeError,
605 "MueLu::Hierarchy::Setup(): fine level (" << startLevel << ") does not exist");
606
607 TEUCHOS_TEST_FOR_EXCEPTION(numDesiredLevels <= 0, Exceptions::RuntimeError,
608 "Constructing non-positive (" << numDesiredLevels << ") number of levels does not make sense.");
609
610 // Check for fine level matrix A
611 TEUCHOS_TEST_FOR_EXCEPTION(!Levels_[startLevel]->IsAvailable("A"), Exceptions::RuntimeError,
612 "MueLu::Hierarchy::Setup(): fine level (" << startLevel << ") has no matrix A! "
613 "Set fine level matrix A using Level.Set()");
614
615 RCP<Operator> A = Levels_[startLevel]->template Get<RCP<Operator> >("A");
616 lib_ = A->getDomainMap()->lib();
617
618 if (IsPrint(Statistics2)) {
619 RCP<Matrix> Amat = rcp_dynamic_cast<Matrix>(A);
620
621 if (!Amat.is_null()) {
622 RCP<ParameterList> params = rcp(new ParameterList());
623 params->set("printLoadBalancingInfo", true);
624 params->set("printCommInfo", true);
625
626 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Amat, "A0", params);
627 } else {
628 GetOStream(Warnings1) << "Fine level operator is not a matrix, statistics are not available" << std::endl;
629 }
630 }
631
632 RCP<const FactoryManagerBase> rcpmanager = rcpFromRef(manager);
633
634 const int lastLevel = startLevel + numDesiredLevels - 1;
635 GetOStream(Runtime0) << "Setup loop: startLevel = " << startLevel << ", lastLevel = " << lastLevel
636 << " (stop if numLevels = " << numDesiredLevels << " or Ac.size() < " << maxCoarseSize_ << ")" << std::endl;
637
638 // Setup multigrid levels
639 int iLevel = 0;
640 if (numDesiredLevels == 1) {
641 iLevel = 0;
642 Setup(startLevel, Teuchos::null, rcpmanager, Teuchos::null); // setup finest==coarsest level (first and last managers are Teuchos::null)
643
644 } else {
645 bool bIsLastLevel = Setup(startLevel, Teuchos::null, rcpmanager, rcpmanager); // setup finest level (level 0) (first manager is Teuchos::null)
646 if (bIsLastLevel == false) {
647 for (iLevel = startLevel + 1; iLevel < lastLevel; iLevel++) {
648 bIsLastLevel = Setup(iLevel, rcpmanager, rcpmanager, rcpmanager); // setup intermediate levels
649 if (bIsLastLevel == true)
650 break;
651 }
652 if (bIsLastLevel == false)
653 Setup(lastLevel, rcpmanager, rcpmanager, Teuchos::null); // setup coarsest level (last manager is Teuchos::null)
654 }
655 }
656
657 // TODO: some check like this should be done at the beginning of the routine
658 TEUCHOS_TEST_FOR_EXCEPTION(iLevel != Levels_.size() - 1, Exceptions::RuntimeError,
659 "MueLu::Hierarchy::Setup(): number of level");
660
661 // TODO: this is not exception safe: manager will still hold default
662 // factories if you exit this function with an exception
663 manager.Clean();
664
666 }
667
668 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
670 if (startLevel < GetNumLevels())
671 GetOStream(Runtime0) << "Clearing old data (if any)" << std::endl;
672
673 for (int iLevel = startLevel; iLevel < GetNumLevels(); iLevel++)
674 Levels_[iLevel]->Clear();
675 }
676
677 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
679 GetOStream(Runtime0) << "Clearing old data (expert)" << std::endl;
680 for (int iLevel = 0; iLevel < GetNumLevels(); iLevel++)
681 Levels_[iLevel]->ExpertClear();
682 }
683
684#if defined(HAVE_MUELU_EXPERIMENTAL) && defined(HAVE_MUELU_ADDITIVE_VARIANT)
685 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
686 ConvergenceStatus Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Iterate(const MultiVector& B, MultiVector& X, ConvData conv,
687 bool InitialGuessIsZero, LO startLevel) {
688 LO nIts = conv.maxIts_;
689 MagnitudeType tol = conv.tol_;
690
691 std::string prefix = this->ShortClassName() + ": ";
692 std::string levelSuffix = " (level=" + toString(startLevel) + ")";
693 std::string levelSuffix1 = " (level=" + toString(startLevel+1) + ")";
694
695 using namespace Teuchos;
696 RCP<Time> CompTime = Teuchos::TimeMonitor::getNewCounter(prefix + "Computational Time (total)");
697 RCP<Time> Concurrent = Teuchos::TimeMonitor::getNewCounter(prefix + "Concurrent portion");
698 RCP<Time> ApplyR = Teuchos::TimeMonitor::getNewCounter(prefix + "R: Computational Time");
699 RCP<Time> ApplyPbar = Teuchos::TimeMonitor::getNewCounter(prefix + "Pbar: Computational Time");
700 RCP<Time> CompFine = Teuchos::TimeMonitor::getNewCounter(prefix + "Fine: Computational Time");
701 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix + "Coarse: Computational Time");
702 RCP<Time> ApplySum = Teuchos::TimeMonitor::getNewCounter(prefix + "Sum: Computational Time");
703 RCP<Time> Synchronize_beginning = Teuchos::TimeMonitor::getNewCounter(prefix + "Synchronize_beginning");
704 RCP<Time> Synchronize_center = Teuchos::TimeMonitor::getNewCounter(prefix + "Synchronize_center");
705 RCP<Time> Synchronize_end = Teuchos::TimeMonitor::getNewCounter(prefix + "Synchronize_end");
706
707 RCP<Level> Fine = Levels_[0];
708 RCP<Level> Coarse;
709
710 RCP<Operator> A = Fine->Get< RCP<Operator> >("A");
711 Teuchos::RCP< const Teuchos::Comm< int > > communicator = A->getDomainMap()->getComm();
712
713
714 //Synchronize_beginning->start();
715 //communicator->barrier();
716 //Synchronize_beginning->stop();
717
718 CompTime->start();
719
720 SC one = STS::one(), zero = STS::zero();
721
722 bool zeroGuess = InitialGuessIsZero;
723
724 // ======= UPFRONT DEFINITION OF COARSE VARIABLES ===========
725
726 //RCP<const Map> origMap;
727 RCP< Operator > P;
728 RCP< Operator > Pbar;
729 RCP< Operator > R;
730 RCP< MultiVector > coarseRhs, coarseX;
731 RCP< Operator > Ac;
732 RCP<SmootherBase> preSmoo_coarse, postSmoo_coarse;
733 bool emptyCoarseSolve = true;
734 RCP<MultiVector> coarseX_prolonged = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(), true);
735
736 RCP<const Import> importer;
737
738 if (Levels_.size() > 1) {
739 Coarse = Levels_[1];
740 if (Coarse->IsAvailable("Importer"))
741 importer = Coarse->Get< RCP<const Import> >("Importer");
742
743 R = Coarse->Get< RCP<Operator> >("R");
744 P = Coarse->Get< RCP<Operator> >("P");
745
746
747 //if(Coarse->IsAvailable("Pbar"))
748 Pbar = Coarse->Get< RCP<Operator> >("Pbar");
749
750 coarseRhs = MultiVectorFactory::Build(R->getRangeMap(), B.getNumVectors(), true);
751
752 Ac = Coarse->Get< RCP< Operator > >("A");
753
754 ApplyR->start();
755 R->apply(B, *coarseRhs, Teuchos::NO_TRANS, one, zero);
756 //P->apply(B, *coarseRhs, Teuchos::TRANS, one, zero);
757 ApplyR->stop();
758
759 if (doPRrebalance_ || importer.is_null()) {
760 coarseX = MultiVectorFactory::Build(coarseRhs->getMap(), X.getNumVectors(), true);
761
762 } else {
763
764 RCP<TimeMonitor> ITime = rcp(new TimeMonitor(*this, prefix + "Solve : import (total)" , Timings0));
765 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : import" + levelSuffix1, Timings0));
766
767 // Import: range map of R --> domain map of rebalanced Ac (before subcomm replacement)
768 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getTargetMap(), coarseRhs->getNumVectors());
769 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
770 coarseRhs.swap(coarseTmp);
771
772 coarseX = MultiVectorFactory::Build(importer->getTargetMap(), X.getNumVectors(), true);
773 }
774
775 if (Coarse->IsAvailable("PreSmoother"))
776 preSmoo_coarse = Coarse->Get< RCP<SmootherBase> >("PreSmoother");
777 if (Coarse->IsAvailable("PostSmoother"))
778 postSmoo_coarse = Coarse->Get< RCP<SmootherBase> >("PostSmoother");
779 }
780
781 // ==========================================================
782
783
784 MagnitudeType prevNorm = STS::magnitude(STS::one()), curNorm = STS::magnitude(STS::one());
785 rate_ = 1.0;
786
787 for (LO i = 1; i <= nIts; i++) {
788#ifdef HAVE_MUELU_DEBUG
789 if (A->getDomainMap()->isCompatible(*(X.getMap())) == false) {
790 std::ostringstream ss;
791 ss << "Level " << startLevel << ": level A's domain map is not compatible with X";
792 throw Exceptions::Incompatible(ss.str());
793 }
794
795 if (A->getRangeMap()->isCompatible(*(B.getMap())) == false) {
796 std::ostringstream ss;
797 ss << "Level " << startLevel << ": level A's range map is not compatible with B";
798 throw Exceptions::Incompatible(ss.str());
799 }
800#endif
801 }
802
803 bool emptyFineSolve = true;
804
805 RCP< MultiVector > fineX;
806 fineX = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(), true);
807
808 //Synchronize_center->start();
809 //communicator->barrier();
810 //Synchronize_center->stop();
811
812 Concurrent->start();
813
814 // NOTE: we need to check using IsAvailable before Get here to avoid building default smoother
815 if (Fine->IsAvailable("PreSmoother")) {
816 RCP<SmootherBase> preSmoo = Fine->Get< RCP<SmootherBase> >("PreSmoother");
817 CompFine->start();
818 preSmoo->Apply(*fineX, B, zeroGuess);
819 CompFine->stop();
820 emptyFineSolve = false;
821 }
822 if (Fine->IsAvailable("PostSmoother")) {
823 RCP<SmootherBase> postSmoo = Fine->Get< RCP<SmootherBase> >("PostSmoother");
824 CompFine->start();
825 postSmoo->Apply(*fineX, B, zeroGuess);
826 CompFine->stop();
827
828 emptyFineSolve = false;
829 }
830 if (emptyFineSolve == true) {
831 GetOStream(Warnings1) << "No fine grid smoother" << std::endl;
832 // Fine grid smoother is identity
833 fineX->update(one, B, zero);
834 }
835
836 if (Levels_.size() > 1) {
837 // NOTE: we need to check using IsAvailable before Get here to avoid building default smoother
838 if (Coarse->IsAvailable("PreSmoother")) {
839 CompCoarse->start();
840 preSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
841 CompCoarse->stop();
842 emptyCoarseSolve = false;
843
844 }
845 if (Coarse->IsAvailable("PostSmoother")) {
846 CompCoarse->start();
847 postSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
848 CompCoarse->stop();
849 emptyCoarseSolve = false;
850
851 }
852 if (emptyCoarseSolve == true) {
853 GetOStream(Warnings1) << "No coarse grid solver" << std::endl;
854 // Coarse operator is identity
855 coarseX->update(one, *coarseRhs, zero);
856 }
857 Concurrent->stop();
858 //Synchronize_end->start();
859 //communicator->barrier();
860 //Synchronize_end->stop();
861
862 if (!doPRrebalance_ && !importer.is_null()) {
863 RCP<TimeMonitor> ITime = rcp(new TimeMonitor(*this, prefix + "Solve : export (total)" , Timings0));
864 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : export" + levelSuffix1, Timings0));
865
866 // Import: range map of rebalanced Ac (before subcomm replacement) --> domain map of P
867 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getSourceMap(), coarseX->getNumVectors());
868 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
869 coarseX.swap(coarseTmp);
870 }
871
872 ApplyPbar->start();
873 Pbar->apply(*coarseX, *coarseX_prolonged, Teuchos::NO_TRANS, one, zero);
874 ApplyPbar->stop();
875 }
876
877 ApplySum->start();
878 X.update(1.0, *fineX, 1.0, *coarseX_prolonged, 0.0);
879 ApplySum->stop();
880
881 CompTime->stop();
882
883 //communicator->barrier();
884
886}
887#else
888 // ---------------------------------------- Iterate -------------------------------------------------------
889 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
891 bool InitialGuessIsZero, LO startLevel) {
892 LO nIts = conv.maxIts_;
893 MagnitudeType tol = conv.tol_;
894
895 // These timers work as follows. "iterateTime" records total time spent in
896 // iterate. "levelTime" records time on a per level basis. The label is
897 // crafted to mimic the per-level messages used in Monitors. Note that a
898 // given level is timed with a TimeMonitor instead of a Monitor or
899 // SubMonitor. This is mainly because I want to time each level
900 // separately, and Monitors/SubMonitors print "(total) xx yy zz" ,
901 // "(sub,total) xx yy zz", respectively, which is subject to
902 // misinterpretation. The per-level TimeMonitors are stopped/started
903 // manually before/after a recursive call to Iterate. A side artifact to
904 // this approach is that the counts for intermediate level timers are twice
905 // the counts for the finest and coarsest levels.
906
907 RCP<Level> Fine = Levels_[startLevel];
908
909 std::string label = FormattingHelper::getColonLabel(Fine->getObjectLabel());
910 std::string prefix = label + this->ShortClassName() + ": ";
911 std::string levelSuffix = " (level=" + toString(startLevel) + ")";
912 std::string levelSuffix1 = " (level=" + toString(startLevel+1) + ")";
913
914 bool useStackedTimer = !Teuchos::TimeMonitor::getStackedTimer().is_null();
915
916 RCP<Monitor> iterateTime;
917 RCP<TimeMonitor> iterateTime1;
918 if (startLevel == 0)
919 iterateTime = rcp(new Monitor(*this, "Solve", label, (nIts == 1) ? None : Runtime0, Timings0));
920 else if (!useStackedTimer)
921 iterateTime1 = rcp(new TimeMonitor(*this, prefix + "Solve (total, level=" + toString(startLevel) + ")", Timings0));
922
923 std::string iterateLevelTimeLabel = prefix + "Solve" + levelSuffix;
924 RCP<TimeMonitor> iterateLevelTime = rcp(new TimeMonitor(*this, iterateLevelTimeLabel, Timings0));
925
926 bool zeroGuess = InitialGuessIsZero;
927
928 RCP<Operator> A = Fine->Get< RCP<Operator> >("A");
929 using namespace Teuchos;
930 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix + "Coarse: Computational Time");
931
932 if (A.is_null()) {
933 // This processor does not have any data for this process on coarser
934 // levels. This can only happen when there are multiple processors and
935 // we use repartitioning.
937 }
938
939 // If we switched the number of vectors, we'd need to reallocate here.
940 // If the number of vectors is unchanged, this is a noop.
941 // NOTE: We need to check against B because the tests in AllocateLevelMultiVectors
942 // will fail on Stokhos Scalar types (due to the so-called 'hidden dimension')
943 const BlockedMultiVector * Bblocked = dynamic_cast<const BlockedMultiVector*>(&B);
944 if(residual_.size() > startLevel &&
945 ( ( Bblocked && !Bblocked->isSameSize(*residual_[startLevel])) ||
946 (!Bblocked && !residual_[startLevel]->isSameSize(B))))
948 AllocateLevelMultiVectors(X.getNumVectors());
949
950 // Print residual information before iterating
951 typedef Teuchos::ScalarTraits<typename STS::magnitudeType> STM;
952 MagnitudeType prevNorm = STM::one();
953 rate_ = 1.0;
954 if (IsCalculationOfResidualRequired(startLevel, conv)) {
955 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, Teuchos::ScalarTraits<LO>::zero(), startLevel, conv, prevNorm);
956 if (convergenceStatus == MueLu::ConvergenceStatus::Converged)
957 return convergenceStatus;
958 }
959
960 SC one = STS::one(), zero = STS::zero();
961 for (LO iteration = 1; iteration <= nIts; iteration++) {
962#ifdef HAVE_MUELU_DEBUG
963#if 0 // TODO fix me
964 if (A->getDomainMap()->isCompatible(*(X.getMap())) == false) {
965 std::ostringstream ss;
966 ss << "Level " << startLevel << ": level A's domain map is not compatible with X";
967 throw Exceptions::Incompatible(ss.str());
968 }
969
970 if (A->getRangeMap()->isCompatible(*(B.getMap())) == false) {
971 std::ostringstream ss;
972 ss << "Level " << startLevel << ": level A's range map is not compatible with B";
973 throw Exceptions::Incompatible(ss.str());
974 }
975#endif
976#endif
977
978 if (startLevel == as<LO>(Levels_.size())-1) {
979 // On the coarsest level, we do either smoothing (if defined) or a direct solve.
980 RCP<TimeMonitor> CLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : coarse" + levelSuffix, Timings0));
981
982 bool emptySolve = true;
983
984 // NOTE: we need to check using IsAvailable before Get here to avoid building default smoother
985 if (Fine->IsAvailable("PreSmoother")) {
986 RCP<SmootherBase> preSmoo = Fine->Get< RCP<SmootherBase> >("PreSmoother");
987 CompCoarse->start();
988 preSmoo->Apply(X, B, zeroGuess);
989 CompCoarse->stop();
990 zeroGuess = false;
991 emptySolve = false;
992 }
993 if (Fine->IsAvailable("PostSmoother")) {
994 RCP<SmootherBase> postSmoo = Fine->Get< RCP<SmootherBase> >("PostSmoother");
995 CompCoarse->start();
996 postSmoo->Apply(X, B, zeroGuess);
997 CompCoarse->stop();
998 emptySolve = false;
999 zeroGuess = false;
1000 }
1001 if (emptySolve == true) {
1002 GetOStream(Warnings1) << "No coarse grid solver" << std::endl;
1003 // Coarse operator is identity
1004 X.update(one, B, zero);
1005 }
1006
1007 } else {
1008 // On intermediate levels, we do cycles
1009 RCP<Level> Coarse = Levels_[startLevel+1];
1010 {
1011 // ============== PRESMOOTHING ==============
1012 RCP<TimeMonitor> STime;
1013 if (!useStackedTimer)
1014 STime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing (total)" , Timings0));
1015 RCP<TimeMonitor> SLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing" + levelSuffix, Timings0));
1016
1017 if (Fine->IsAvailable("PreSmoother")) {
1018 RCP<SmootherBase> preSmoo = Fine->Get< RCP<SmootherBase> >("PreSmoother");
1019 preSmoo->Apply(X, B, zeroGuess);
1020 zeroGuess = false;
1021 }
1022 }
1023
1024 RCP<MultiVector> residual;
1025 {
1026 RCP<TimeMonitor> ATime;
1027 if (!useStackedTimer)
1028 ATime = rcp(new TimeMonitor(*this, prefix + "Solve : residual calculation (total)" , Timings0));
1029 RCP<TimeMonitor> ALevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : residual calculation" + levelSuffix, Timings0));
1030 if (zeroGuess) {
1031 // If there's a pre-smoother, then zeroGuess is false. If there isn't and people still have zeroGuess set,
1032 // then then X still has who knows what, so we need to zero it before we go to the coarse grid.
1033 X.putScalar(zero);
1034 }
1035
1036 Utilities::Residual(*A, X, B,*residual_[startLevel]);
1037 residual = residual_[startLevel];
1038 }
1039
1040 RCP<Operator> P = Coarse->Get< RCP<Operator> >("P");
1041 if (Coarse->IsAvailable("Pbar"))
1042 P = Coarse->Get< RCP<Operator> >("Pbar");
1043
1044 RCP<MultiVector> coarseRhs, coarseX;
1045 // const bool initializeWithZeros = true;
1046 {
1047 // ============== RESTRICTION ==============
1048 RCP<TimeMonitor> RTime;
1049 if (!useStackedTimer)
1050 RTime = rcp(new TimeMonitor(*this, prefix + "Solve : restriction (total)" , Timings0));
1051 RCP<TimeMonitor> RLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : restriction" + levelSuffix, Timings0));
1052 coarseRhs = coarseRhs_[startLevel];
1053
1054 if (implicitTranspose_) {
1055 P->apply(*residual, *coarseRhs, Teuchos::TRANS, one, zero);
1056
1057 } else {
1058 RCP<Operator> R = Coarse->Get< RCP<Operator> >("R");
1059 R->apply(*residual, *coarseRhs, Teuchos::NO_TRANS, one, zero);
1060 }
1061 }
1062
1063 RCP<const Import> importer;
1064 if (Coarse->IsAvailable("Importer"))
1065 importer = Coarse->Get< RCP<const Import> >("Importer");
1066
1067 coarseX = coarseX_[startLevel];
1068 if (!doPRrebalance_ && !importer.is_null()) {
1069 RCP<TimeMonitor> ITime;
1070 if (!useStackedTimer)
1071 ITime = rcp(new TimeMonitor(*this, prefix + "Solve : import (total)" , Timings0));
1072 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : import" + levelSuffix1, Timings0));
1073
1074 // Import: range map of R --> domain map of rebalanced Ac (before subcomm replacement)
1075 RCP<MultiVector> coarseTmp = coarseImport_[startLevel];
1076 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
1077 coarseRhs.swap(coarseTmp);
1078 }
1079
1080 RCP<Operator> Ac = Coarse->Get< RCP<Operator> >("A");
1081 if (!Ac.is_null()) {
1082 RCP<const Map> origXMap = coarseX->getMap();
1083 RCP<const Map> origRhsMap = coarseRhs->getMap();
1084
1085 // Replace maps with maps with a subcommunicator
1086 coarseRhs->replaceMap(Ac->getRangeMap());
1087 coarseX ->replaceMap(Ac->getDomainMap());
1088
1089 {
1090 iterateLevelTime = Teuchos::null; // stop timing this level
1091
1092 Iterate(*coarseRhs, *coarseX, 1, true, startLevel+1);
1093 // ^^ zero initial guess
1094 if (Cycle_ == WCYCLE && WCycleStartLevel_ >= startLevel)
1095 Iterate(*coarseRhs, *coarseX, 1, false, startLevel+1);
1096 // ^^ nonzero initial guess
1097
1098 iterateLevelTime = rcp(new TimeMonitor(*this, iterateLevelTimeLabel)); // restart timing this level
1099 }
1100 coarseX->replaceMap(origXMap);
1101 coarseRhs->replaceMap(origRhsMap);
1102 }
1103
1104 if (!doPRrebalance_ && !importer.is_null()) {
1105 RCP<TimeMonitor> ITime;
1106 if (!useStackedTimer)
1107 ITime = rcp(new TimeMonitor(*this, prefix + "Solve : export (total)" , Timings0));
1108 RCP<TimeMonitor> ILevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : export" + levelSuffix1, Timings0));
1109
1110 // Import: range map of rebalanced Ac (before subcomm replacement) --> domain map of P
1111 RCP<MultiVector> coarseTmp = coarseExport_[startLevel];
1112 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
1113 coarseX.swap(coarseTmp);
1114 }
1115
1116 {
1117 // ============== PROLONGATION ==============
1118 RCP<TimeMonitor> PTime;
1119 if (!useStackedTimer)
1120 PTime = rcp(new TimeMonitor(*this, prefix + "Solve : prolongation (total)" , Timings0));
1121 RCP<TimeMonitor> PLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : prolongation" + levelSuffix, Timings0));
1122 // Update X += P * coarseX
1123 // Note that due to what may be round-off error accumulation, use of the fused kernel
1124 // P->apply(*coarseX, X, Teuchos::NO_TRANS, one, one);
1125 // can in some cases result in slightly higher iteration counts.
1127 P->apply(*coarseX, X, Teuchos::NO_TRANS, scalingFactor_, one);
1128 } else {
1129 RCP<MultiVector> correction = correction_[startLevel];
1130 P->apply(*coarseX, *correction, Teuchos::NO_TRANS, one, zero);
1131 X.update(scalingFactor_, *correction, one);
1132 }
1133 }
1134
1135 {
1136 // ============== POSTSMOOTHING ==============
1137 RCP<TimeMonitor> STime;
1138 if (!useStackedTimer)
1139 STime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing (total)" , Timings0));
1140 RCP<TimeMonitor> SLevelTime = rcp(new TimeMonitor(*this, prefix + "Solve : smoothing" + levelSuffix, Timings0));
1141
1142 if (Fine->IsAvailable("PostSmoother")) {
1143 RCP<SmootherBase> postSmoo = Fine->Get< RCP<SmootherBase> >("PostSmoother");
1144 postSmoo->Apply(X, B, false);
1145 }
1146 }
1147 }
1148 zeroGuess = false;
1149
1150
1151 if (IsCalculationOfResidualRequired(startLevel, conv)) {
1152 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, iteration, startLevel, conv, prevNorm);
1153 if (convergenceStatus == MueLu::ConvergenceStatus::Converged)
1154 return convergenceStatus;
1155 }
1156 }
1158 }
1159#endif
1160
1161
1162 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1163 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(const LO& start, const LO& end, const std::string &suffix) {
1164 LO startLevel = (start != -1 ? start : 0);
1165 LO endLevel = (end != -1 ? end : Levels_.size()-1);
1166
1167 TEUCHOS_TEST_FOR_EXCEPTION(startLevel > endLevel, Exceptions::RuntimeError,
1168 "MueLu::Hierarchy::Write : startLevel must be <= endLevel");
1169
1170 TEUCHOS_TEST_FOR_EXCEPTION(startLevel < 0 || endLevel >= Levels_.size(), Exceptions::RuntimeError,
1171 "MueLu::Hierarchy::Write bad start or end level");
1172
1173 for (LO i = startLevel; i < endLevel + 1; i++) {
1174 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Levels_[i]-> template Get< RCP< Operator> >("A")), P, R;
1175 if (i > 0) {
1176 P = rcp_dynamic_cast<Matrix>(Levels_[i]-> template Get< RCP< Operator> >("P"));
1177 if (!implicitTranspose_)
1178 R = rcp_dynamic_cast<Matrix>(Levels_[i]-> template Get< RCP< Operator> >("R"));
1179 }
1180
1181 if (!A.is_null()) Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write("A_" + toString(i) + suffix + ".m", *A);
1182 if (!P.is_null()) {
1183 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write("P_" + toString(i) + suffix + ".m", *P);
1184 }
1185 if (!R.is_null()) {
1186 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write("R_" + toString(i) + suffix + ".m", *R);
1187 }
1188 }
1189 }
1190
1191 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1192 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Keep(const std::string & ename, const FactoryBase* factory) {
1193 for (Array<RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1194 (*it)->Keep(ename, factory);
1195 }
1196
1197 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1198 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Delete(const std::string& ename, const FactoryBase* factory) {
1199 for (Array<RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1200 (*it)->Delete(ename, factory);
1201 }
1202
1203 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1205 for (Array<RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1206 (*it)->AddKeepFlag(ename, factory, keep);
1207 }
1208
1209 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1211 for (Array<RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1212 (*it)->RemoveKeepFlag(ename, factory, keep);
1213 }
1214
1215 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1217 if ( description_.empty() )
1218 {
1219 std::ostringstream out;
1220 out << BaseClass::description();
1221 out << "{#levels = " << GetGlobalNumLevels() << ", complexity = " << GetOperatorComplexity() << "}";
1222 description_ = out.str();
1223 }
1224 return description_;
1225 }
1226
1227 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1228 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel tVerbLevel) const {
1229 describe(out, toMueLuVerbLevel(tVerbLevel));
1230 }
1231
1232 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1233 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::describe(Teuchos::FancyOStream& out, const VerbLevel verbLevel) const {
1234 RCP<Operator> A0 = Levels_[0]->template Get<RCP<Operator> >("A");
1235 RCP<const Teuchos::Comm<int> > comm = A0->getDomainMap()->getComm();
1236
1237 int numLevels = GetNumLevels();
1238 RCP<Operator> Ac = Levels_[numLevels-1]->template Get<RCP<Operator> >("A");
1239 if (Ac.is_null()) {
1240 // It may happen that we do repartition on the last level, but the matrix
1241 // is small enough to satisfy "max coarse size" requirement. Then, even
1242 // though we have the level, the matrix would be null on all but one processors
1243 numLevels--;
1244 }
1245 int root = comm->getRank();
1246
1247#ifdef HAVE_MPI
1248 int smartData = numLevels*comm->getSize() + comm->getRank(), maxSmartData;
1249 reduceAll(*comm, Teuchos::REDUCE_MAX, smartData, Teuchos::ptr(&maxSmartData));
1250 root = maxSmartData % comm->getSize();
1251#endif
1252
1253 // Compute smoother complexity, if needed
1254 double smoother_comp =-1.0;
1255 if (verbLevel & (Statistics0 | Test))
1256 smoother_comp = GetSmootherComplexity();
1257
1258 std::string outstr;
1259 if (comm->getRank() == root && verbLevel & (Statistics0 | Test)) {
1260 std::vector<Xpetra::global_size_t> nnzPerLevel;
1261 std::vector<Xpetra::global_size_t> rowsPerLevel;
1262 std::vector<int> numProcsPerLevel;
1263 bool someOpsNotMatrices = false;
1264 const Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1265 for (int i = 0; i < numLevels; i++) {
1266 TEUCHOS_TEST_FOR_EXCEPTION(!(Levels_[i]->IsAvailable("A")) , Exceptions::RuntimeError,
1267 "Operator A is not available on level " << i);
1268
1269 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator> >("A");
1270 TEUCHOS_TEST_FOR_EXCEPTION(A.is_null(), Exceptions::RuntimeError,
1271 "Operator A on level " << i << " is null.");
1272
1273 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
1274 if (Am.is_null()) {
1275 someOpsNotMatrices = true;
1276 nnzPerLevel .push_back(INVALID);
1277 rowsPerLevel .push_back(A->getDomainMap()->getGlobalNumElements());
1278 numProcsPerLevel.push_back(A->getDomainMap()->getComm()->getSize());
1279 } else {
1280 LO storageblocksize=Am->GetStorageBlockSize();
1281 Xpetra::global_size_t nnz = Am->getGlobalNumEntries()*storageblocksize*storageblocksize;
1282 nnzPerLevel .push_back(nnz);
1283 rowsPerLevel .push_back(Am->getGlobalNumRows()*storageblocksize);
1284 numProcsPerLevel.push_back(Am->getRowMap()->getComm()->getSize());
1285 }
1286 }
1287 if (someOpsNotMatrices)
1288 GetOStream(Warnings0) << "Some level operators are not matrices, statistics calculation are incomplete" << std::endl;
1289
1290 {
1291 std::string label = Levels_[0]->getObjectLabel();
1292 std::ostringstream oss;
1293 oss << std::setfill(' ');
1294 oss << "\n--------------------------------------------------------------------------------\n";
1295 oss << "--- Multigrid Summary " << std::setw(28) << std::left << label << "---\n";
1296 oss << "--------------------------------------------------------------------------------" << std::endl;
1297 if (verbLevel & Parameters1)
1298 oss << "Scalar = " << Teuchos::ScalarTraits<Scalar>::name() << std::endl;
1299 oss << "Number of levels = " << numLevels << std::endl;
1300 oss << "Operator complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1301 if (!someOpsNotMatrices)
1302 oss << GetOperatorComplexity() << std::endl;
1303 else
1304 oss << "not available (Some operators in hierarchy are not matrices.)" << std::endl;
1305
1306 if(smoother_comp!=-1.0) {
1307 oss << "Smoother complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed)
1308 << smoother_comp << std::endl;
1309 }
1310
1311 switch (Cycle_) {
1312 case VCYCLE:
1313 oss << "Cycle type = V" << std::endl;
1314 break;
1315 case WCYCLE:
1316 oss << "Cycle type = W" << std::endl;
1317 if (WCycleStartLevel_ > 0)
1318 oss << "Cycle start level = " << WCycleStartLevel_ << std::endl;
1319 break;
1320 default:
1321 break;
1322 };
1323 oss << std::endl;
1324
1325 Xpetra::global_size_t tt = rowsPerLevel[0];
1326 int rowspacer = 2; while (tt != 0) { tt /= 10; rowspacer++; }
1327 for (size_t i = 0; i < nnzPerLevel.size(); ++i) {
1328 tt = nnzPerLevel[i];
1329 if (tt != INVALID)
1330 break;
1331 tt = 100; // This will get used if all levels are operators.
1332 }
1333 int nnzspacer = 2; while (tt != 0) { tt /= 10; nnzspacer++; }
1334 tt = numProcsPerLevel[0];
1335 int npspacer = 2; while (tt != 0) { tt /= 10; npspacer++; }
1336 oss << "level " << std::setw(rowspacer) << " rows " << std::setw(nnzspacer) << " nnz " << " nnz/row" << std::setw(npspacer) << " c ratio" << " procs" << std::endl;
1337 for (size_t i = 0; i < nnzPerLevel.size(); ++i) {
1338 oss << " " << i << " ";
1339 oss << std::setw(rowspacer) << rowsPerLevel[i];
1340 if (nnzPerLevel[i] != INVALID) {
1341 oss << std::setw(nnzspacer) << nnzPerLevel[i];
1342 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1343 oss << std::setw(9) << as<double>(nnzPerLevel[i]) / rowsPerLevel[i];
1344 } else {
1345 oss << std::setw(nnzspacer) << "Operator";
1346 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1347 oss << std::setw(9) << " ";
1348 }
1349 if (i) oss << std::setw(9) << as<double>(rowsPerLevel[i-1])/rowsPerLevel[i];
1350 else oss << std::setw(9) << " ";
1351 oss << " " << std::setw(npspacer) << numProcsPerLevel[i] << std::endl;
1352 }
1353 oss << std::endl;
1354 for (int i = 0; i < GetNumLevels(); ++i) {
1355 RCP<SmootherBase> preSmoo, postSmoo;
1356 if (Levels_[i]->IsAvailable("PreSmoother"))
1357 preSmoo = Levels_[i]->template Get< RCP<SmootherBase> >("PreSmoother");
1358 if (Levels_[i]->IsAvailable("PostSmoother"))
1359 postSmoo = Levels_[i]->template Get< RCP<SmootherBase> >("PostSmoother");
1360
1361 if (preSmoo != null && preSmoo == postSmoo)
1362 oss << "Smoother (level " << i << ") both : " << preSmoo->description() << std::endl;
1363 else {
1364 oss << "Smoother (level " << i << ") pre : "
1365 << (preSmoo != null ? preSmoo->description() : "no smoother") << std::endl;
1366 oss << "Smoother (level " << i << ") post : "
1367 << (postSmoo != null ? postSmoo->description() : "no smoother") << std::endl;
1368 }
1369
1370 oss << std::endl;
1371 }
1372
1373 outstr = oss.str();
1374 }
1375 }
1376
1377#ifdef HAVE_MPI
1378 RCP<const Teuchos::MpiComm<int> > mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int> >(comm);
1379 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
1380
1381 int strLength = outstr.size();
1382 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
1383 if (comm->getRank() != root)
1384 outstr.resize(strLength);
1385 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
1386#endif
1387
1388 out << outstr;
1389 }
1390
1391 // NOTE: at some point this should be replaced by a friend operator <<
1392 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1393 void Hierarchy<Scalar, LocalOrdinal, GlobalOrdinal, Node>::print(std::ostream& out, const VerbLevel verbLevel) const {
1394 Teuchos::OSTab tab2(out);
1395 for (int i = 0; i < GetNumLevels(); ++i)
1396 Levels_[i]->print(out, verbLevel);
1397 }
1398
1399 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1403
1404 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1406 if (GetProcRankVerbose() != 0)
1407 return;
1408#if defined(HAVE_MUELU_BOOST) && defined(HAVE_MUELU_BOOST_FOR_REAL) && defined(BOOST_VERSION) && (BOOST_VERSION >= 104400)
1409
1410 BoostGraph graph;
1411
1412 BoostProperties dp;
1413 dp.property("label", boost::get(boost::vertex_name, graph));
1414 dp.property("id", boost::get(boost::vertex_index, graph));
1415 dp.property("label", boost::get(boost::edge_name, graph));
1416 dp.property("color", boost::get(boost::edge_color, graph));
1417
1418 // create local maps
1419 std::map<const FactoryBase*, BoostVertex> vindices;
1420 typedef std::map<std::pair<BoostVertex,BoostVertex>, std::string> emap; emap edges;
1421
1422 static int call_id=0;
1423
1424 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator> >("A");
1425 int rank = A->getDomainMap()->getComm()->getRank();
1426
1427 // printf("[%d] CMS: ----------------------\n",rank);
1428 for (int i = currLevel; i <= currLevel+1 && i < GetNumLevels(); i++) {
1429 edges.clear();
1430 Levels_[i]->UpdateGraph(vindices, edges, dp, graph);
1431
1432 for (emap::const_iterator eit = edges.begin(); eit != edges.end(); eit++) {
1433 std::pair<BoostEdge, bool> boost_edge = boost::add_edge(eit->first.first, eit->first.second, graph);
1434 // printf("[%d] CMS: Hierarchy, adding edge (%d->%d) %d\n",rank,(int)eit->first.first,(int)eit->first.second,(int)boost_edge.second);
1435 // Because xdot.py views 'Graph' as a keyword
1436 if(eit->second==std::string("Graph")) boost::put("label", dp, boost_edge.first, std::string("Graph_"));
1437 else boost::put("label", dp, boost_edge.first, eit->second);
1438 if (i == currLevel)
1439 boost::put("color", dp, boost_edge.first, std::string("red"));
1440 else
1441 boost::put("color", dp, boost_edge.first, std::string("blue"));
1442 }
1443 }
1444
1445 std::ofstream out(dumpFile_.c_str()+std::string("_")+std::to_string(currLevel)+std::string("_")+std::to_string(call_id)+std::string("_")+ std::to_string(rank) + std::string(".dot"));
1446 boost::write_graphviz_dp(out, graph, dp, std::string("id"));
1447 out.close();
1448 call_id++;
1449#else
1450 GetOStream(Errors) << "Dependency graph output requires boost and MueLu_ENABLE_Boost_for_real" << std::endl;
1451#endif
1452 }
1453
1454 // Enforce that coordinate vector's map is consistent with that of A
1455 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1457 RCP<Operator> Ao = level.Get<RCP<Operator> >("A");
1458 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Ao);
1459 if (A.is_null()) {
1460 GetOStream(Runtime1) << "Hierarchy::ReplaceCoordinateMap: operator is not a matrix, skipping..." << std::endl;
1461 return;
1462 }
1463 if(Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A) != Teuchos::null) {
1464 GetOStream(Runtime1) << "Hierarchy::ReplaceCoordinateMap: operator is a BlockedCrsMatrix, skipping..." << std::endl;
1465 return;
1466 }
1467
1468 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO> xdMV;
1469
1470 RCP<xdMV> coords = level.Get<RCP<xdMV> >("Coordinates");
1471
1472 if (A->getRowMap()->isSameAs(*(coords->getMap()))) {
1473 GetOStream(Runtime1) << "Hierarchy::ReplaceCoordinateMap: matrix and coordinates maps are same, skipping..." << std::endl;
1474 return;
1475 }
1476
1477 if (A->IsView("stridedMaps") && rcp_dynamic_cast<const StridedMap>(A->getRowMap("stridedMaps")) != Teuchos::null) {
1478 RCP<const StridedMap> stridedRowMap = rcp_dynamic_cast<const StridedMap>(A->getRowMap("stridedMaps"));
1479
1480 // It is better to through an exceptions if maps may be inconsistent, than to ignore it and experience unfathomable breakdowns
1481 TEUCHOS_TEST_FOR_EXCEPTION(stridedRowMap->getStridedBlockId() != -1 || stridedRowMap->getOffset() != 0,
1482 Exceptions::RuntimeError, "Hierarchy::ReplaceCoordinateMap: nontrivial maps (block id = " << stridedRowMap->getStridedBlockId()
1483 << ", offset = " << stridedRowMap->getOffset() << ")");
1484 }
1485
1486 GetOStream(Runtime1) << "Replacing coordinate map" << std::endl;
1487 TEUCHOS_TEST_FOR_EXCEPTION(A->GetFixedBlockSize() % A->GetStorageBlockSize() != 0, Exceptions::RuntimeError, "Hierarchy::ReplaceCoordinateMap: Storage block size does not evenly divide fixed block size");
1488
1489 size_t blkSize = A->GetFixedBlockSize() / A->GetStorageBlockSize();
1490
1491 RCP<const Map> nodeMap = A->getRowMap();
1492 if (blkSize > 1) {
1493 // Create a nodal map, as coordinates have not been expanded to a DOF map yet.
1494 RCP<const Map> dofMap = A->getRowMap();
1495 GO indexBase = dofMap->getIndexBase();
1496 size_t numLocalDOFs = dofMap->getLocalNumElements();
1497 TEUCHOS_TEST_FOR_EXCEPTION(numLocalDOFs % blkSize, Exceptions::RuntimeError,
1498 "Hierarchy::ReplaceCoordinateMap: block size (" << blkSize << ") is incompatible with the number of local dofs in a row map (" << numLocalDOFs);
1499 ArrayView<const GO> GIDs = dofMap->getLocalElementList();
1500
1501 Array<GO> nodeGIDs(numLocalDOFs/blkSize);
1502 for (size_t i = 0; i < numLocalDOFs; i += blkSize)
1503 nodeGIDs[i/blkSize] = (GIDs[i] - indexBase)/blkSize + indexBase;
1504
1505 Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1506 nodeMap = MapFactory::Build(dofMap->lib(), INVALID, nodeGIDs(), indexBase, dofMap->getComm());
1507 } else {
1508 // blkSize == 1
1509 // Check whether the length of vectors fits to the size of A
1510 // If yes, make sure that the maps are matching
1511 // If no, throw a warning but do not touch the Coordinates
1512 if(coords->getLocalLength() != A->getRowMap()->getLocalNumElements()) {
1513 GetOStream(Warnings) << "Coordinate vector does not match row map of matrix A!" << std::endl;
1514 return;
1515 }
1516 }
1517
1518 Array<ArrayView<const typename Teuchos::ScalarTraits<Scalar>::coordinateType> > coordDataView;
1519 std::vector<ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType> > coordData;
1520 for (size_t i = 0; i < coords->getNumVectors(); i++) {
1521 coordData.push_back(coords->getData(i));
1522 coordDataView.push_back(coordData[i]());
1523 }
1524
1525 RCP<xdMV> newCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType,LO,GO,NO>::Build(nodeMap, coordDataView(), coords->getNumVectors());
1526 level.Set("Coordinates", newCoords);
1527 }
1528
1529 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1531 int N = Levels_.size();
1532 if( ( (sizeOfAllocatedLevelMultiVectors_ == numvecs && residual_.size() == N) || numvecs<=0 ) && !forceMapCheck) return;
1533
1534 // If, somehow, we changed the number of levels, delete everything first
1535 if(residual_.size() != N) {
1537
1538 residual_.resize(N);
1539 coarseRhs_.resize(N);
1540 coarseX_.resize(N);
1541 coarseImport_.resize(N);
1542 coarseExport_.resize(N);
1543 correction_.resize(N);
1544 }
1545
1546 for(int i=0; i<N; i++) {
1547 RCP<Operator> A = Levels_[i]->template Get< RCP<Operator> >("A");
1548 if(!A.is_null()) {
1549 // This dance is because we allow A to have a BlockedMap and X/B to have (compatible) non-blocked map
1550 RCP<const BlockedCrsMatrix> A_as_blocked = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(A);
1551 RCP<const Map> Arm = A->getRangeMap();
1552 RCP<const Map> Adm = A->getDomainMap();
1553 if(!A_as_blocked.is_null()) {
1554 Adm = A_as_blocked->getFullDomainMap();
1555 }
1556
1557 if (residual_[i].is_null() || !residual_[i]->getMap()->isSameAs(*Arm))
1558 // This is zero'd by default since it is filled via an operator apply
1559 residual_[i] = MultiVectorFactory::Build(Arm, numvecs, true);
1560 if (correction_[i].is_null() || !correction_[i]->getMap()->isSameAs(*Adm))
1561 correction_[i] = MultiVectorFactory::Build(Adm, numvecs, false);
1562 }
1563
1564 if(i+1<N) {
1565 // This is zero'd by default since it is filled via an operator apply
1566 if(implicitTranspose_) {
1567 RCP<Operator> P = Levels_[i+1]->template Get< RCP<Operator> >("P");
1568 if(!P.is_null()) {
1569 RCP<const Map> map = P->getDomainMap();
1570 if (coarseRhs_[i].is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1571 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs, true);
1572 }
1573 } else {
1574 RCP<Operator> R = Levels_[i+1]->template Get< RCP<Operator> >("R");
1575 if (!R.is_null()) {
1576 RCP<const Map> map = R->getRangeMap();
1577 if (coarseRhs_[i].is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1578 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs, true);
1579 }
1580 }
1581
1582
1583 RCP<const Import> importer;
1584 if(Levels_[i+1]->IsAvailable("Importer"))
1585 importer = Levels_[i+1]->template Get< RCP<const Import> >("Importer");
1586 if (doPRrebalance_ || importer.is_null()) {
1587 RCP<const Map> map = coarseRhs_[i]->getMap();
1588 if (coarseX_[i].is_null() || !coarseX_[i]->getMap()->isSameAs(*map))
1589 coarseX_[i] = MultiVectorFactory::Build(map, numvecs, true);
1590 } else {
1591 RCP<const Map> map;
1592 map = importer->getTargetMap();
1593 if (coarseImport_[i].is_null() || !coarseImport_[i]->getMap()->isSameAs(*map)) {
1594 coarseImport_[i] = MultiVectorFactory::Build(map, numvecs,false);
1595 coarseX_[i] = MultiVectorFactory::Build(map,numvecs,false);
1596 }
1597 map = importer->getSourceMap();
1598 if (coarseExport_[i].is_null() || !coarseExport_[i]->getMap()->isSameAs(*map))
1599 coarseExport_[i] = MultiVectorFactory::Build(map, numvecs,false);
1600 }
1601 }
1602 }
1604 }
1605
1606
1607template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1610 residual_.resize(0);
1611 coarseRhs_.resize(0);
1612 coarseX_.resize(0);
1613 coarseImport_.resize(0);
1614 coarseExport_.resize(0);
1615 correction_.resize(0);
1617}
1618
1619
1620template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1622 const LO startLevel, const ConvData& conv) const
1623{
1624 return (startLevel == 0 && !isPreconditioner_ && (IsPrint(Statistics1) || conv.tol_ > 0));
1625}
1626
1627
1628template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1630 const Teuchos::Array<MagnitudeType>& residualNorm, const MagnitudeType convergenceTolerance) const
1631{
1633
1634 if (convergenceTolerance > Teuchos::ScalarTraits<MagnitudeType>::zero())
1635 {
1636 bool passed = true;
1637 for (LO k = 0; k < residualNorm.size(); k++)
1638 if (residualNorm[k] >= convergenceTolerance)
1639 passed = false;
1640
1641 if (passed)
1642 convergenceStatus = ConvergenceStatus::Converged;
1643 else
1644 convergenceStatus = ConvergenceStatus::Unconverged;
1645 }
1646
1647 return convergenceStatus;
1648}
1649
1650
1651template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1653 const LO iteration, const Teuchos::Array<MagnitudeType>& residualNorm) const
1654{
1655 GetOStream(Statistics1) << "iter: "
1656 << std::setiosflags(std::ios::left)
1657 << std::setprecision(3) << std::setw(4) << iteration
1658 << " residual = "
1659 << std::setprecision(10) << residualNorm
1660 << std::endl;
1661}
1662
1663template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1665 const Operator& A, const MultiVector& X, const MultiVector& B, const LO iteration,
1666 const LO startLevel, const ConvData& conv, MagnitudeType& previousResidualNorm)
1667{
1668 Teuchos::Array<MagnitudeType> residualNorm;
1669 residualNorm = Utilities::ResidualNorm(A, X, B, *residual_[startLevel]);
1670
1671 const MagnitudeType currentResidualNorm = residualNorm[0];
1672 rate_ = currentResidualNorm / previousResidualNorm;
1673 previousResidualNorm = currentResidualNorm;
1674
1675 if (IsPrint(Statistics1))
1676 PrintResidualHistory(iteration, residualNorm);
1677
1678 return IsConverged(residualNorm, conv.tol_);
1679}
1680
1681
1682} //namespace MueLu
1683
1684#endif // MUELU_HIERARCHY_DEF_HPP
virtual std::string ShortClassName() const
Return the class name of the object, without template parameters and without namespace.
virtual std::string description() const
Return a simple one-line description of this object.
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
Class that provides default factories within Needs class.
void AddLevel(const RCP< Level > &level)
Add a level at the end of the hierarchy.
double GetSmootherComplexity() const
void Write(const LO &start=-1, const LO &end=-1, const std::string &suffix="")
Print matrices in the multigrid hierarchy to file.
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
Array< RCP< MultiVector > > residual_
void CheckLevel(Level &level, int levelID)
Helper function.
int WCycleStartLevel_
Level at which to start W-cycle.
std::string description() const
Return a simple one-line description of this object.
std::string description_
cache description to avoid recreating in each call to description() - use ResetDescription() to force...
void IsPreconditioner(const bool flag)
static CycleType GetDefaultCycle()
Array< RCP< Level > > Levels_
Container for Level objects.
static bool GetDefaultPRrebalance()
Array< RCP< MultiVector > > coarseRhs_
Array< RCP< MultiVector > > coarseExport_
bool Setup(int coarseLevelID, const RCP< const FactoryManagerBase > fineLevelManager, const RCP< const FactoryManagerBase > coarseLevelManager, const RCP< const FactoryManagerBase > nextLevelManager=Teuchos::null)
Multi-level setup phase: build a new level of the hierarchy.
void ResetDescription()
force recreation of cached description_ next time description() is called:
STS::magnitudeType MagnitudeType
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
bool isPreconditioner_
Hierarchy may be used in a standalone mode, or as a preconditioner.
MueLu::Level Level
ConvergenceStatus IsConverged(const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const
Decide if the multigrid iteration is converged.
ConvergenceStatus Iterate(const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0)
Apply the multigrid preconditioner.
void DumpCurrentGraph(int level) const
int sizeOfAllocatedLevelMultiVectors_
Caching (Multi)Vectors used in Hierarchy::Iterate().
static bool GetDefaultImplicitTranspose()
void SetMatvecParams(RCP< ParameterList > matvecParams)
static Xpetra::global_size_t GetDefaultMaxCoarseSize()
Xpetra::UnderlyingLib lib_
Epetra/Tpetra mode.
void Clear(int startLevel=0)
Clear impermanent data from previous setup.
bool IsCalculationOfResidualRequired(const LO startLevel, const ConvData &conv) const
Decide if the residual needs to be computed.
CycleType Cycle_
V- or W-cycle.
ConvergenceStatus ComputeResidualAndPrintHistory(const Operator &A, const MultiVector &X, const MultiVector &B, const LO iteration, const LO startLevel, const ConvData &conv, MagnitudeType &previousResidualNorm)
Compute the residual norm and print it depending on the verbosity level.
double GetOperatorComplexity() const
MagnitudeType rate_
Convergece rate.
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
Array< RCP< MultiVector > > coarseX_
void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck=false)
void print(std::ostream &out=std::cout, const VerbLevel verbLevel=(MueLu::Parameters|MueLu::Statistics0)) const
Hierarchy::print is local hierarchy function, thus the statistics can be different from global ones.
double scalingFactor_
Scaling factor to be applied to coarse grid correction.
void Delete(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Delete(ename, factory) for each level of the Hierarchy.
void Keep(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Keep(ename, factory) for each level of the Hierarchy.
Xpetra::UnderlyingLib lib()
static bool GetDefaultFuseProlongationAndUpdate()
Array< RCP< MultiVector > > coarseImport_
bool isDumpingEnabled_
Graph dumping.
Xpetra::global_size_t maxCoarseSize_
MueLu::FactoryBase FactoryBase
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
Array< RCP< MultiVector > > correction_
void AddNewLevel()
Add a new level at the end of the hierarchy.
void RemoveKeepFlag(const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All)
Call Level::RemoveKeepFlag for each level of the Hierarchy.
Array< RCP< const FactoryManagerBase > > levelManagers_
Level managers used during the Setup.
void ReplaceCoordinateMap(Level &level)
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void SetComm(RCP< const Teuchos::Comm< int > > const &comm)
void setlib(Xpetra::UnderlyingLib lib2)
RCP< Level > & GetPreviousLevel()
Previous level.
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
RCP< const Teuchos::Comm< int > > GetComm() const
int GetLevelID() const
Return level number.
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())
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.
Xpetra::UnderlyingLib lib()
Timer to be used in non-factories.
static const NoFactory * get()
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
An exception safe way to call the method 'Level::SetFactoryManager()'.
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
VerbLevel GetVerbLevel() const
Get the verbosity level.
int SetProcRankVerbose(int procRank) const
Set proc rank used for printing.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
int GetProcRankVerbose() const
Get proc rank used for printing. Do not use this information for any other purpose....
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line).
@ Statistics2
Print even more statistics.
@ Warnings
Print all warning messages.
@ Statistics1
Print more statistics.
@ Timings0
High level timing information (use Teuchos::TimeMonitor::summarize() to print).
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose).
@ Warnings1
Additional warnings.
@ Statistics0
Print statistics that do not involve significant additional computation.
@ Parameters1
Print class parameters (more parameters, more verbose).
short KeepType
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
VerbLevel toMueLuVerbLevel(const Teuchos::EVerbosityLevel verbLevel)
Translate Teuchos verbosity level to MueLu verbosity level.
static std::string getColonLabel(const std::string &label)
Helper function for object label.
Data struct for defining stopping criteria of multigrid iteration.