MueLu  Version of the Day
MueLu_RebalanceBlockAcFactory_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_REBALANCEBLOCKACFACTORY_DEF_HPP_
48 #define MUELU_REBALANCEBLOCKACFACTORY_DEF_HPP_
49 
50 #include <Xpetra_Matrix.hpp>
51 #include <Xpetra_CrsMatrix.hpp>
52 #include <Xpetra_CrsMatrixWrap.hpp>
53 #include <Xpetra_MatrixFactory.hpp>
54 #include <Xpetra_MapExtractor.hpp>
55 #include <Xpetra_MapExtractorFactory.hpp>
56 #include <Xpetra_StridedMap.hpp>
57 #include <Xpetra_StridedMapFactory.hpp>
58 #include <Xpetra_BlockedCrsMatrix.hpp>
59 
60 #include <Xpetra_VectorFactory.hpp>
61 
63 
65 #include "MueLu_HierarchyUtils.hpp"
66 #include "MueLu_MasterList.hpp"
67 #include "MueLu_Monitor.hpp"
68 #include "MueLu_PerfUtils.hpp"
69 #include "MueLu_RAPFactory.hpp"
70 
71 namespace MueLu {
72 
73  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75  RCP<ParameterList> validParamList = rcp(new ParameterList());
76 
77 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
78  SET_VALID_ENTRY("repartition: use subcommunicators");
79 #undef SET_VALID_ENTRY
80 
81  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A for rebalancing");
82  validParamList->set<RCP<const FactoryBase> >("Importer", Teuchos::null, "Generating factory of the matrix Importer for rebalancing");
83  validParamList->set<RCP<const FactoryBase> >("SubImporters", Teuchos::null, "Generating factory of the matrix sub-Importers for rebalancing");
84 
85  return validParamList;
86  }
87 
88  template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
90  FactManager_.push_back(FactManager);
91  }
92 
93  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95  Input(coarseLevel, "A");
96 
97  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
98  for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
99  SetFactoryManager fineSFM (rcpFromRef(fineLevel), *it);
100  SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
101 
102  coarseLevel.DeclareInput("Importer",(*it)->GetFactory("Importer").get(), this);
103  }
104 
105  // Use the non-manager path if the maps / importers are generated in one place
106  if(FactManager_.size() == 0) {
107  Input(coarseLevel,"SubImporters");
108  }
109 
110  }
111 
112  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
114  FactoryMonitor m(*this, "Computing blocked Ac", coarseLevel);
115 
116  RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
117 
118  RCP<Matrix> originalAc = Get< RCP<Matrix> >(coarseLevel, "A");
119 
120  RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > bA = Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(originalAc);
121  TEUCHOS_TEST_FOR_EXCEPTION(bA==Teuchos::null, Exceptions::BadCast, "MueLu::RebalanceBlockAcFactory::Build: input matrix A is not of type BlockedCrsMatrix! error.");
122  TEUCHOS_TEST_FOR_EXCEPTION(bA->Rows() != bA->Cols(), Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: Blocked operator has " << bA->Rows() << " and " << bA->Cols() << ". We only support square matrices (with same number of blocks and columns).");
123 
124  // Variables to set up map extractors for blocked operators
125  std::vector<GO> fullRangeMapVector;
126  std::vector<GO> fullDomainMapVector;
127  std::vector<RCP<const Map> > subBlockARangeMaps;
128  std::vector<RCP<const Map> > subBlockADomainMaps;
129  subBlockARangeMaps.reserve(bA->Rows());
130  subBlockADomainMaps.reserve(bA->Cols());
131 
132  // store map extractors
133  RCP<const MapExtractor> rangeMapExtractor = bA->getRangeMapExtractor();
134  RCP<const MapExtractor> domainMapExtractor = bA->getDomainMapExtractor();
135 
136  // check if GIDs for full maps have to be sorted:
137  // For the Thyra mode ordering they do not have to be sorted since the GIDs are
138  // numbered as 0...n1,0...,n2 (starting with zero for each subblock). The MapExtractor
139  // generates unique GIDs during the construction.
140  // For Xpetra style, the GIDs have to be reordered. Such that one obtains a ordered
141  // list of GIDs in an increasing ordering. In Xpetra, the GIDs are all unique through
142  // out all submaps.
143  bool bThyraRangeGIDs = rangeMapExtractor->getThyraMode();
144  bool bThyraDomainGIDs = domainMapExtractor->getThyraMode();
145 
146  // vector containing rebalanced blocks (for final output)
147  std::vector<RCP<Matrix> > subBlockRebA =
148  std::vector<RCP<Matrix> >(bA->Cols() * bA->Rows(), Teuchos::null);
149 
150  // vector with Import objects from the different
151  // RepartitionFactory instances
152  std::vector<RCP<const Import> > importers = std::vector<RCP<const Import> >(bA->Rows(), Teuchos::null);
153  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
154  size_t idx = 0;
155  for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
156  SetFactoryManager fineSFM (rcpFromRef(fineLevel), *it);
157  SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
158 
159  RCP<const Import> rebalanceImporter = coarseLevel.Get<RCP<const Import> >("Importer", (*it)->GetFactory("Importer").get());
160  importers[idx] = rebalanceImporter;
161  idx++;
162  }
163  if(FactManager_.size() == 0) {
164  importers = Get<std::vector<RCP<const Import> > >(coarseLevel,"SubImporters");
165  }
166 
167 
168  // restrict communicator?
169  bool bRestrictComm = false;
170  const ParameterList& pL = GetParameterList();
171  if (pL.get<bool>("repartition: use subcommunicators") == true)
172  bRestrictComm = true;
173 
174  RCP<ParameterList> XpetraList = Teuchos::rcp(new ParameterList());
175  if(bRestrictComm)
176  XpetraList->set("Restrict Communicator",true);
177  else
178  XpetraList->set("Restrict Communicator",false);
179 
180  // communicator for final (rebalanced) operator.
181  // If communicator is not restricted it should be identical to the communicator in bA
182  RCP<const Teuchos::Comm<int> > rebalancedComm = Teuchos::null;
183 
184  // loop through all blocks and rebalance blocks
185  // Note: so far we do not support rebalancing of nested operators
186  // TODO add a check for this
187  for(size_t i=0; i<bA->Rows(); i++) {
188  for(size_t j=0; j<bA->Cols(); j++) {
189  // extract matrix block
190  RCP<Matrix> Aij = bA->getMatrix(i, j);
191 
192  std::stringstream ss; ss << "Rebalancing matrix block A(" << i << "," << j << ")";
193  SubFactoryMonitor subM(*this, ss.str(), coarseLevel);
194 
195  RCP<Matrix> rebAij = Teuchos::null;
196  // General rebalancing
197  if( importers[i] != Teuchos::null &&
198  importers[j] != Teuchos::null &&
199  Aij != Teuchos::null) {
200  RCP<const Map> targetRangeMap = importers[i]->getTargetMap();
201  RCP<const Map> targetDomainMap = importers[j]->getTargetMap();
202 
203  // Copy the block Aij
204  // TAW: Do we need a copy or can we do in-place rebalancing?
205  // If we do in-place rebalancing the original distribution is lost
206  // We don't really need it any more, though.
207  //RCP<Matrix> cAij = MatrixFactory::BuildCopy(Aij);
208  RCP<Matrix> cAij = Aij; // do not copy the matrix data (just an rcp pointer)
209 
210  // create a new importer for column map needed for rebalanced Aij
211  Teuchos::RCP<const Import> rebAijImport = ImportFactory::Build(importers[j]->getTargetMap(),cAij->getColMap());
212  TEUCHOS_TEST_FOR_EXCEPTION(rebAijImport.is_null() == true,Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: Importer associated with block " << j << " is null.");
213 
214  Teuchos::RCP<const CrsMatrixWrap> cAwij = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(cAij);
215  TEUCHOS_TEST_FOR_EXCEPTION(cAwij.is_null() == true,Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: Block (" << i << "," << j << ") is not of type CrsMatrix. We cannot rebalanced (nested) operators.");
216  Teuchos::RCP<CrsMatrix> cAmij = cAwij->getCrsMatrix();
217 
218  // change domain map to rebalanced domain map (in-place). Update the importer to represent the column map
219  //cAmij->replaceDomainMapAndImporter(importers[j]->getTargetMap(),rebAijImport);
220 
221  // rebalance rows of matrix block. Don't change the domain map (-> Teuchos::null)
222  // NOTE: If the communicator is restricted away, Build returns Teuchos::null.
223  rebAij = MatrixFactory::Build(cAij, *(importers[i]), *(importers[j]), targetDomainMap, targetRangeMap, XpetraList);
224  } // rebalance matrix block A(i,i)
225  else {
226  rebAij = Aij; // no rebalancing or empty block!
227  }
228 
229  // store new block in output
230  subBlockRebA[i*bA->Cols() + j] = rebAij;
231 
232  if (!rebAij.is_null()) {
233  // store communicator
234  if(rebalancedComm.is_null()) rebalancedComm = rebAij->getRowMap()->getComm();
235 
236  // printout rebalancing information
237  RCP<ParameterList> params = rcp(new ParameterList());
238  params->set("printLoadBalancingInfo", true);
239  std::stringstream ss2; ss2 << "A(" << i << "," << j << ") rebalanced:";
240  GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*rebAij, ss2.str(), params);
241  }
242  } // loop over columns j
243 
244  // fix striding information of diagonal blocks
245  // Note: we do not care about the off-diagonal blocks. We just make sure, that the
246  // diagonal blocks have the corresponding striding information from the map extractors
247  // Note: the diagonal block never should be zero.
248  // TODO what if a diagonal block is Teuchos::null?
249  if ( subBlockRebA[i*bA->Cols() + i].is_null() == false ) {
250  RCP<Matrix> rebAii = subBlockRebA[i*bA->Cols() + i];
251  Teuchos::RCP<const StridedMap> orig_stridedRgMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getMap(i,rangeMapExtractor->getThyraMode()));
252  Teuchos::RCP<const Map> stridedRgMap = Teuchos::null;
253  if(orig_stridedRgMap != Teuchos::null) {
254  std::vector<size_t> stridingData = orig_stridedRgMap->getStridingData();
255  Teuchos::ArrayView< const GlobalOrdinal > nodeRangeMapii = rebAii->getRangeMap()->getNodeElementList();
256  stridedRgMap = StridedMapFactory::Build(
257  bA->getRangeMap()->lib(),
258  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
259  nodeRangeMapii,
260  rebAii->getRangeMap()->getIndexBase(),
261  stridingData,
262  rebalancedComm,
263  orig_stridedRgMap->getStridedBlockId(),
264  orig_stridedRgMap->getOffset());
265  }
266  Teuchos::RCP<const StridedMap> orig_stridedDoMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getMap(i,domainMapExtractor->getThyraMode()));
267  Teuchos::RCP<const Map> stridedDoMap = Teuchos::null;
268  if(orig_stridedDoMap != Teuchos::null) {
269  std::vector<size_t> stridingData = orig_stridedDoMap->getStridingData();
270  Teuchos::ArrayView< const GlobalOrdinal > nodeDomainMapii = rebAii->getDomainMap()->getNodeElementList();
271  stridedDoMap = StridedMapFactory::Build(
272  bA->getDomainMap()->lib(),
273  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
274  nodeDomainMapii,
275  rebAii->getDomainMap()->getIndexBase(),
276  stridingData,
277  rebalancedComm,
278  orig_stridedDoMap->getStridedBlockId(),
279  orig_stridedDoMap->getOffset());
280  }
281 
282  if(bRestrictComm) {
283  stridedRgMap->removeEmptyProcesses();
284  stridedDoMap->removeEmptyProcesses();
285  }
286 
287  TEUCHOS_TEST_FOR_EXCEPTION(stridedRgMap == Teuchos::null,Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: failed to generate striding information. error.");
288  TEUCHOS_TEST_FOR_EXCEPTION(stridedDoMap == Teuchos::null,Exceptions::RuntimeError, "MueLu::RebalanceBlockAcFactory::Build: failed to generate striding information. error.");
289 
290  // replace stridedMaps view in diagonal sub block
291  if(rebAii->IsView("stridedMaps")) rebAii->RemoveView("stridedMaps");
292  rebAii->CreateView("stridedMaps", stridedRgMap, stridedDoMap);
293  // collect Xpetra-based global row ids for map extractors
294  subBlockARangeMaps.push_back(rebAii->getRowMap("stridedMaps"));
295  Teuchos::ArrayView< const GlobalOrdinal > nodeRangeMap = rebAii->getRangeMap()->getNodeElementList();
296  // append the GIDs in the end. Do not sort if we have Thyra style GIDs
297  fullRangeMapVector.insert(fullRangeMapVector.end(), nodeRangeMap.begin(), nodeRangeMap.end());
298  if(bThyraRangeGIDs == false)
299  sort(fullRangeMapVector.begin(), fullRangeMapVector.end());
300 
301  subBlockADomainMaps.push_back(rebAii->getColMap("stridedMaps"));
302  Teuchos::ArrayView< const GlobalOrdinal > nodeDomainMap = rebAii->getDomainMap()->getNodeElementList();
303  // append the GIDs in the end. Do not sort if we have Thyra style GIDs
304  fullDomainMapVector.insert(fullDomainMapVector.end(), nodeDomainMap.begin(), nodeDomainMap.end());
305  if(bThyraDomainGIDs == false)
306  sort(fullDomainMapVector.begin(), fullDomainMapVector.end());
307  } // end if rebAii != Teuchos::null
308  } // loop over rows i
309 
310  // all sub blocks are rebalanced (if available)
311 
312  // Short cut if this processor is not in the list of active processors
313  if (rebalancedComm == Teuchos::null) {
314  GetOStream(Debug,-1) << "RebalanceBlockedAc: deactivate proc " << originalAc->getRowMap()->getComm()->getRank() << std::endl;
315  // TAW: it is important that we create a dummy object of type BlockedCrsMatrix (even if we set it to Teuchos::null)
316  Teuchos::RCP<BlockedCrsMatrix> reb_bA = Teuchos::null;
317  coarseLevel.Set("A", Teuchos::rcp_dynamic_cast<Matrix>(reb_bA), this);
318  return;
319  }
320 
321  // now, subBlockRebA contains all rebalanced matrix blocks
322  // extract map index base from maps of blocked A
323  GO rangeIndexBase = bA->getRangeMap()->getIndexBase();
324  GO domainIndexBase = bA->getDomainMap()->getIndexBase();
325 
326  Teuchos::ArrayView<GO> fullRangeMapGIDs(fullRangeMapVector.size() ? &fullRangeMapVector[0] : 0,fullRangeMapVector.size());
327  Teuchos::RCP<const StridedMap> stridedRgFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getFullMap());
328  Teuchos::RCP<const Map > fullRangeMap = Teuchos::null;
329  if(stridedRgFullMap != Teuchos::null) {
330  std::vector<size_t> stridedData = stridedRgFullMap->getStridingData();
331  fullRangeMap =
332  StridedMapFactory::Build(
333  bA->getRangeMap()->lib(),
334  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
335  fullRangeMapGIDs,
336  rangeIndexBase,
337  stridedData,
338  rebalancedComm,
339  stridedRgFullMap->getStridedBlockId(),
340  stridedRgFullMap->getOffset());
341  } else {
342  fullRangeMap =
343  MapFactory::Build(
344  bA->getRangeMap()->lib(),
345  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
346  fullRangeMapGIDs,
347  rangeIndexBase,
348  rebalancedComm);
349  }
350  Teuchos::ArrayView<GO> fullDomainMapGIDs(fullDomainMapVector.size() ? &fullDomainMapVector[0] : 0,fullDomainMapVector.size());
351 
352  Teuchos::RCP<const StridedMap> stridedDoFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getFullMap());
353  Teuchos::RCP<const Map > fullDomainMap = Teuchos::null;
354  if(stridedDoFullMap != Teuchos::null) {
355  TEUCHOS_TEST_FOR_EXCEPTION(stridedDoFullMap==Teuchos::null, Exceptions::BadCast, "MueLu::RebalanceBlockedAc::Build: full map in domain map extractor has no striding information! error.");
356  std::vector<size_t> stridedData2 = stridedDoFullMap->getStridingData();
357  fullDomainMap =
358  StridedMapFactory::Build(
359  bA->getDomainMap()->lib(),
360  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
361  fullDomainMapGIDs,
362  domainIndexBase,
363  stridedData2,
364  rebalancedComm,
365  stridedDoFullMap->getStridedBlockId(),
366  stridedDoFullMap->getOffset());
367  } else {
368 
369  fullDomainMap =
370  MapFactory::Build(
371  bA->getDomainMap()->lib(),
372  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
373  fullDomainMapGIDs,
374  domainIndexBase,
375  rebalancedComm);
376  }
377 
378  if(bRestrictComm) {
379  fullRangeMap->removeEmptyProcesses();
380  fullDomainMap->removeEmptyProcesses();
381  }
382 
383  // build map extractors
384  RCP<const MapExtractor> rebRangeMapExtractor = MapExtractorFactory::Build(fullRangeMap, subBlockARangeMaps, bThyraRangeGIDs);
385  RCP<const MapExtractor> rebDomainMapExtractor = MapExtractorFactory::Build(fullDomainMap, subBlockADomainMaps, bThyraDomainGIDs);
386 
387  TEUCHOS_TEST_FOR_EXCEPTION(rangeMapExtractor->NumMaps() != rebRangeMapExtractor->NumMaps(), Exceptions::RuntimeError,
388  "MueLu::RebalanceBlockedAc::Build: Rebalanced RangeMapExtractor has " << rebRangeMapExtractor->NumMaps()
389  << " sub maps. Original RangeMapExtractor has " << rangeMapExtractor->NumMaps() << ". They must match!");
390  TEUCHOS_TEST_FOR_EXCEPTION(domainMapExtractor->NumMaps() != rebDomainMapExtractor->NumMaps(), Exceptions::RuntimeError,
391  "MueLu::RebalanceBlockedAc::Build: Rebalanced DomainMapExtractor has " << rebDomainMapExtractor->NumMaps()
392  << " sub maps. Original DomainMapExtractor has " << domainMapExtractor->NumMaps() << ". They must match!");
393 
394  Teuchos::RCP<BlockedCrsMatrix> reb_bA = Teuchos::rcp(new BlockedCrsMatrix(rebRangeMapExtractor,rebDomainMapExtractor,10));
395  for(size_t i=0; i<bA->Rows(); i++) {
396  for(size_t j=0; j<bA->Cols(); j++) {
397  reb_bA->setMatrix(i,j,subBlockRebA[i*bA->Cols() + j]);
398  }
399  }
400 
401  reb_bA->fillComplete();
402 
403  coarseLevel.Set("A", Teuchos::rcp_dynamic_cast<Matrix>(reb_bA), this);
404  // rebalance additional data:
405  // be aware, that we just call the rebalance factories without switching to local
406  // factory managers, i.e. the rebalance factories have to be defined with the appropriate
407  // factories by the user!
408  if (rebalanceFacts_.begin() != rebalanceFacts_.end()) {
409  SubFactoryMonitor m2(*this, "Rebalance additional data", coarseLevel);
410 
411  // call Build of all user-given transfer factories
412  for (std::vector<RCP<const FactoryBase> >::const_iterator it2 = rebalanceFacts_.begin(); it2 != rebalanceFacts_.end(); ++it2) {
413  GetOStream(Runtime0) << "RebalanceBlockedAc: call rebalance factory " << (*it2).get() << ": " << (*it2)->description() << std::endl;
414  (*it2)->CallBuild(coarseLevel);
415  }
416  }
417  } //Build()
418 
419  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
421  rebalanceFacts_.push_back(factory);
422  } //AddRebalanceFactory()
423 
424 } //namespace MueLu
425 
426 #endif /* MUELU_REBALANCEBLOCKACFACTORY_DEF_HPP_ */
#define SET_VALID_ENTRY(name)
Exception indicating invalid cast attempted.
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.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
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())
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager)
Add a factory manager.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void AddRebalanceFactory(const RCP< const FactoryBase > &factory)
Add rebalancing factory in the end of list of rebalancing factories in RebalanceAcFactory.
An exception safe way to call the method 'Level::SetFactoryManager()'.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Debug
Print additional debugging information.
@ Runtime0
One-liner description of what is happening.
@ Statistics0
Print statistics that do not involve significant additional computation.