MueLu  Version of the Day
MueLu_RebalanceBlockInterpolationFactory_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 #ifndef MUELU_REBALANCEBLOCKINTERPOLATIONFACTORY_DEF_HPP_
47 #define MUELU_REBALANCEBLOCKINTERPOLATIONFACTORY_DEF_HPP_
48 
49 #ifdef HAVE_MUELU_EXPERIMENTAL
50 
51 #include <Teuchos_Tuple.hpp>
52 
53 #include "Xpetra_MultiVector.hpp"
54 #include "Xpetra_MultiVectorFactory.hpp"
55 #include "Xpetra_Vector.hpp"
56 #include "Xpetra_VectorFactory.hpp"
57 #include <Xpetra_Matrix.hpp>
58 #include <Xpetra_BlockedCrsMatrix.hpp>
59 #include <Xpetra_MapFactory.hpp>
60 #include <Xpetra_MapExtractor.hpp>
61 #include <Xpetra_MapExtractorFactory.hpp>
62 #include <Xpetra_MatrixFactory.hpp>
63 #include <Xpetra_Import.hpp>
64 #include <Xpetra_ImportFactory.hpp>
65 
67 
69 #include "MueLu_HierarchyUtils.hpp"
70 #include "MueLu_Level.hpp"
71 #include "MueLu_Monitor.hpp"
72 #include "MueLu_PerfUtils.hpp"
73 
74 namespace MueLu {
75 
76 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78  RCP<ParameterList> validParamList = rcp(new ParameterList());
79 
80  validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
81  validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Factory for generating the non-rebalanced coarse level A. We need this to make sure the non-rebalanced coarse A is calculated first before rebalancing takes place.");
82 
83  validParamList->set< RCP<const FactoryBase> >("Coordinates", Teuchos::null, "Factory for generating the non-rebalanced Coordinates.");
84 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
85  // SET_VALID_ENTRY("repartition: use subcommunicators");
86 #undef SET_VALID_ENTRY
87 
88  // TODO validation: "P" parameter valid only for type="Interpolation" and "R" valid only for type="Restriction". Like so:
89  // if (paramList.isEntry("type") && paramList.get("type) == "Interpolation) {
90  // validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Factory of the prolongation operator that need to be rebalanced (only used if type=Interpolation)");
91 
92  return validParamList;
93 }
94 
95 template <class Scalar,class LocalOrdinal, class GlobalOrdinal, class Node>
97  FactManager_.push_back(FactManager);
98 }
99 
100 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
102  Input(coarseLevel, "P");
103  Input(coarseLevel, "A"); // we request the non-rebalanced coarse level A since we have to make sure it is calculated before rebalancing starts!
104 
105  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
106  for(it = FactManager_.begin(); it!=FactManager_.end(); ++it) {
107  SetFactoryManager fineSFM (rcpFromRef(fineLevel), *it);
108  SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
109 
110  // Request Importer and Coordinates (if defined in xml file)
111  // Note, that we have to use the Level::DeclareInput routine in order to use the FactoryManager *it (rather than the main factory manager)
112  coarseLevel.DeclareInput("Importer",(*it)->GetFactory("Importer").get(), this);
113  if((*it)->hasFactory("Coordinates") == true)
114  coarseLevel.DeclareInput("Coordinates",(*it)->GetFactory("Coordinates").get(), this);
115  }
116 }
117 
118 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
120  FactoryMonitor m(*this, "Build", coarseLevel);
121 
122  //RCP<Teuchos::FancyOStream> out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
123 
124  Teuchos::RCP<Matrix> nonrebCoarseA = Get< RCP<Matrix> >(coarseLevel, "A");
125 
126  Teuchos::RCP<Matrix> originalTransferOp = Teuchos::null;
127  originalTransferOp = Get< RCP<Matrix> >(coarseLevel, "P");
128 
129  RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > bOriginalTransferOp =
130  Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > (originalTransferOp);
131  TEUCHOS_TEST_FOR_EXCEPTION(bOriginalTransferOp==Teuchos::null, Exceptions::BadCast, "MueLu::RebalanceBlockTransferFactory::Build: input matrix P or R is not of type BlockedCrsMatrix! error.");
132 
133  RCP<const MapExtractor> rangeMapExtractor = bOriginalTransferOp->getRangeMapExtractor();
134  RCP<const MapExtractor> domainMapExtractor = bOriginalTransferOp->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  // declare variables for maps of blocked rebalanced prolongation operator
147  std::vector<GO> fullRangeMapVector; // contains all range GIDs on current processor
148  std::vector<GO> fullDomainMapVector; // contains all domain GIDs on current processor
149  std::vector<RCP<const Map> > subBlockPRangeMaps;
150  std::vector<RCP<const Map> > subBlockPDomainMaps;
151  subBlockPRangeMaps.reserve(bOriginalTransferOp->Rows()); // reserve size for block P operators
152  subBlockPDomainMaps.reserve(bOriginalTransferOp->Cols()); // reserve size for block P operators
153 
154  std::vector<Teuchos::RCP<Matrix> > subBlockRebP;
155  subBlockRebP.reserve(bOriginalTransferOp->Rows());
156 
157  int curBlockId = 0;
158  Teuchos::RCP<const Import> rebalanceImporter = Teuchos::null;
159  std::vector<Teuchos::RCP<const FactoryManagerBase> >::const_iterator it;
160  for (it = FactManager_.begin(); it != FactManager_.end(); ++it) {
161  // begin SubFactoryManager environment
162  SetFactoryManager fineSFM (rcpFromRef(fineLevel), *it);
163  SetFactoryManager coarseSFM(rcpFromRef(coarseLevel), *it);
164 
165  // TAW: use the Level::Get routine in order to access the data declared in (*it) factory manager (rather than the main factory manager)
166  rebalanceImporter = coarseLevel.Get<Teuchos::RCP<const Import> >("Importer", (*it)->GetFactory("Importer").get());
167 
168  // extract diagonal matrix block
169  Teuchos::RCP<Matrix> Pii = bOriginalTransferOp->getMatrix(curBlockId, curBlockId);
170  Teuchos::RCP<CrsMatrixWrap> Pwii = Teuchos::rcp_dynamic_cast<CrsMatrixWrap>(Pii);
171  TEUCHOS_TEST_FOR_EXCEPTION(Pwii == Teuchos::null,Xpetra::Exceptions::BadCast, "MueLu::RebalanceBlockTransferFactory::Build: block " << curBlockId << " is not of type CrsMatrixWrap. We need an underlying CsrMatrix to replace domain map and importer!");
172 
173  MUELU_TEST_FOR_EXCEPTION(bThyraRangeGIDs == true && Pii->getRowMap()->getMinAllGlobalIndex() != 0,
175  "MueLu::RebalanceBlockInterpolationFactory::Build: inconsistent Thyra GIDs. Thyra global ids for block range " << curBlockId << " start with " << Pii->getRowMap()->getMinAllGlobalIndex() << " but should start with 0");
176  MUELU_TEST_FOR_EXCEPTION(bThyraDomainGIDs == true && Pii->getColMap()->getMinAllGlobalIndex() != 0,
178  "MueLu::RebalanceBlockInterpolationFactory::Build: inconsistent Thyra GIDs. Thyra global ids for block domain " << curBlockId << " start with " << Pii->getColMap()->getMinAllGlobalIndex() << " but should start with 0");
179 
180  // rebalance P11
181  if(rebalanceImporter != Teuchos::null) {
182  std::stringstream ss; ss << "Rebalancing prolongator block P(" << curBlockId << "," << curBlockId << ")";
183  SubFactoryMonitor m1(*this, ss.str(), coarseLevel);
184 
185  // P is the transfer operator from the coarse grid to the fine grid.
186  // P must transfer the data from the newly reordered coarse A to the (unchanged) fine A.
187  // This means that the domain map (coarse) of P must be changed according to the new partition. The range map (fine) is kept unchanged.
188  //
189  // The domain map of P must match the range map of R.
190  // See also note below about domain/range map of R and its implications for P.
191  //
192  // To change the domain map of P, P needs to be fillCompleted again with the new domain map.
193  // To achieve this, P is copied into a new matrix that is not fill-completed.
194  // The doImport() operation is just used here to make a copy of P: the importer is trivial and there is no data movement involved.
195  // The reordering actually happens during the fillComplete() with domainMap == rebalanceImporter->getTargetMap().
196  RCP<const Import> newImporter;
197  {
198  SubFactoryMonitor subM(*this, "Rebalancing prolongator -- fast map replacement", coarseLevel);
199  newImporter = ImportFactory::Build(rebalanceImporter->getTargetMap(), Pii->getColMap());
200  Pwii->getCrsMatrix()->replaceDomainMapAndImporter(rebalanceImporter->getTargetMap(), newImporter);
201  }
202 
203  RCP<ParameterList> params = rcp(new ParameterList());
204  params->set("printLoadBalancingInfo", true);
205  std::stringstream ss2; ss2 << "P(" << curBlockId << "," << curBlockId << ") rebalanced:";
206  GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*Pii, ss2.str(), params);
207 
208  // store rebalanced P block
209  subBlockRebP.push_back(Pii);
210 
211  // rebalance coordinates
212  // TAW: Note, that each sub-block manager overwrites the Coordinates. So far we only support one set of Coordinates
213  // for a multiphysics problem (i.e., we only support volume coupled problems with the same mesh)
214 
215  if((*it)->hasFactory("Coordinates") == true && coarseLevel.IsAvailable("Coordinates",(*it)->GetFactory("Coordinates").get()) == true) {
216  typedef Xpetra::MultiVector<double, LO, GO, NO> xdMV;
217  RCP<xdMV> coords = coarseLevel.Get< RCP<xdMV> >("Coordinates",(*it)->GetFactory("Coordinates").get());
218 
219  // This line must be after the Get call
220  SubFactoryMonitor subM(*this, "Rebalancing coordinates", coarseLevel);
221 
222  LO nodeNumElts = coords->getMap()->getNodeNumElements();
223 
224  // If a process has no matrix rows, then we can't calculate blocksize using the formula below.
225  LO myBlkSize = 0, blkSize = 0;
226 
227  if (nodeNumElts > 0) {
228  MUELU_TEST_FOR_EXCEPTION(rebalanceImporter->getSourceMap()->getNodeNumElements() % nodeNumElts != 0,
230  "MueLu::RebalanceBlockInterpolationFactory::Build: block size. " << rebalanceImporter->getSourceMap()->getNodeNumElements() << " not divisable by " << nodeNumElts);
231  myBlkSize = rebalanceImporter->getSourceMap()->getNodeNumElements() / nodeNumElts;
232  }
233 
234  MueLu_maxAll(coords->getMap()->getComm(), myBlkSize, blkSize);
235 
236  RCP<const Import> coordImporter = Teuchos::null;
237  if (blkSize == 1) {
238  coordImporter = rebalanceImporter;
239 
240  } else {
241  // NOTE: there is an implicit assumption here: we assume that dof any node are enumerated consequently
242  // Proper fix would require using decomposition similar to how we construct importer in the
243  // RepartitionFactory
244  RCP<const Map> origMap = coords->getMap();
245  GO indexBase = origMap->getIndexBase();
246 
247  ArrayView<const GO> OEntries = rebalanceImporter->getTargetMap()->getNodeElementList();
248  LO numEntries = OEntries.size()/blkSize;
249  ArrayRCP<GO> Entries(numEntries);
250  for (LO i = 0; i < numEntries; i++)
251  Entries[i] = (OEntries[i*blkSize]-indexBase)/blkSize + indexBase;
252 
253  RCP<const Map> targetMap = MapFactory::Build(origMap->lib(), origMap->getGlobalNumElements(), Entries(), indexBase, origMap->getComm());
254  coordImporter = ImportFactory::Build(origMap, targetMap);
255  }
256 
257  RCP<xdMV> permutedCoords = Xpetra::MultiVectorFactory<double,LO,GO,NO>::Build(coordImporter->getTargetMap(), coords->getNumVectors());
258  permutedCoords->doImport(*coords, *coordImporter, Xpetra::INSERT);
259 
260  const ParameterList& pL = GetParameterList();
261  if (pL.isParameter("repartition: use subcommunicators") == true && pL.get<bool>("repartition: use subcommunicators") == true)
262  permutedCoords->replaceMap(permutedCoords->getMap()->removeEmptyProcesses());
263 
264  Set(coarseLevel, "Coordinates", permutedCoords);
265  }
266  } // end rebalance P(1,1)
267  else {
268  RCP<ParameterList> params = rcp(new ParameterList());
269  params->set("printLoadBalancingInfo", true);
270  std::stringstream ss; ss << "P(" << curBlockId << "," << curBlockId << ") not rebalanced:";
271  GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*Pii, ss.str(), params);
272  // store rebalanced P block
273  subBlockRebP.push_back(Pii);
274 
275  // Store Coordinates on coarse level (generated by this)
276  // TAW: Note, that each sub-block manager overwrites the Coordinates. So far we only support one set of Coordinates
277  // for a multiphysics problem (i.e., we only support volume coupled problems with the same mesh)
278  if((*it)->hasFactory("Coordinates") == true && coarseLevel.IsAvailable("Coordinates",(*it)->GetFactory("Coordinates").get()) == true) {
279  typedef Xpetra::MultiVector<double, LO, GO, NO> xdMV;
280  coarseLevel.Set("Coordinates", coarseLevel.Get< RCP<xdMV> >("Coordinates",(*it)->GetFactory("Coordinates").get()),this);
281  }
282  }
283 
284  // fix striding information for rebalanced diagonal block Pii
285  //RCP<const Xpetra::MapExtractor<Scalar, LocalOrdinal, GlobalOrdinal, Node> > rgPMapExtractor = bOriginalTransferOp->getRangeMapExtractor(); // original map extractor
286  Teuchos::RCP<const StridedMap> orig_stridedRgMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangeMapExtractor->getMap(Teuchos::as<size_t>(curBlockId),rangeMapExtractor->getThyraMode()));
287  Teuchos::RCP<const Map> stridedRgMap = Teuchos::null;
288  if(orig_stridedRgMap != Teuchos::null) {
289  std::vector<size_t> stridingData = orig_stridedRgMap->getStridingData();
290  Teuchos::ArrayView< const GlobalOrdinal > nodeRangeMapii = Pii->getRangeMap()->getNodeElementList();
291  stridedRgMap = StridedMapFactory::Build(
292  originalTransferOp->getRangeMap()->lib(),
293  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
294  nodeRangeMapii,
295  Pii->getRangeMap()->getIndexBase(),
296  stridingData,
297  originalTransferOp->getRangeMap()->getComm(),
298  orig_stridedRgMap->getStridedBlockId(),
299  orig_stridedRgMap->getOffset());
300  } else stridedRgMap = Pii->getRangeMap();
301 
302  //RCP<const Xpetra::MapExtractor<Scalar, LocalOrdinal, GlobalOrdinal, Node> > doPMapExtractor = bOriginalTransferOp->getDomainMapExtractor(); // original map extractor
303  Teuchos::RCP<const StridedMap> orig_stridedDoMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainMapExtractor->getMap(Teuchos::as<size_t>(curBlockId),domainMapExtractor->getThyraMode()));
304 
305  Teuchos::RCP<const Map> stridedDoMap = Teuchos::null;
306  if(orig_stridedDoMap != Teuchos::null) {
307  std::vector<size_t> stridingData = orig_stridedDoMap->getStridingData();
308  Teuchos::ArrayView< const GlobalOrdinal > nodeDomainMapii = Pii->getDomainMap()->getNodeElementList();
309  stridedDoMap = StridedMapFactory::Build(originalTransferOp->getDomainMap()->lib(),
310  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
311  nodeDomainMapii,
312  Pii->getDomainMap()->getIndexBase(),
313  stridingData,
314  originalTransferOp->getDomainMap()->getComm(),
315  orig_stridedDoMap->getStridedBlockId(),
316  orig_stridedDoMap->getOffset());
317  } else stridedDoMap = Pii->getDomainMap();
318 
319  TEUCHOS_TEST_FOR_EXCEPTION(stridedRgMap == Teuchos::null,Exceptions::RuntimeError, "MueLu::RebalanceBlockInterpolationFactory::Build: failed to generate striding information. error.");
320  TEUCHOS_TEST_FOR_EXCEPTION(stridedDoMap == Teuchos::null,Exceptions::RuntimeError, "MueLu::RebalanceBlockInterpolationFactory::Build: failed to generate striding information. error.");
321 
322  // replace stridedMaps view in diagonal sub block
323  if(Pii->IsView("stridedMaps")) Pii->RemoveView("stridedMaps");
324  Pii->CreateView("stridedMaps", stridedRgMap, stridedDoMap);
325 
326  // append strided row map (= range map) to list of range maps.
327  Teuchos::RCP<const Map> rangeMapii = Pii->getRowMap("stridedMaps"); // strided range map
328  subBlockPRangeMaps.push_back(rangeMapii);
329  Teuchos::ArrayView< const GlobalOrdinal > nodeRangeMapii = Pii->getRangeMap()->getNodeElementList();
330  // append the GIDs in the end. Do not sort if we have Thyra style GIDs
331  fullRangeMapVector.insert(fullRangeMapVector.end(), nodeRangeMapii.begin(), nodeRangeMapii.end());
332  if(bThyraRangeGIDs == false)
333  sort(fullRangeMapVector.begin(), fullRangeMapVector.end());
334 
335  // append strided col map (= domain map) to list of range maps.
336  Teuchos::RCP<const Map> domainMapii = Pii->getColMap("stridedMaps"); // strided domain map
337  subBlockPDomainMaps.push_back(domainMapii);
338  Teuchos::ArrayView< const GlobalOrdinal > nodeDomainMapii = Pii->getDomainMap()->getNodeElementList();
339  // append the GIDs in the end. Do not sort if we have Thyra style GIDs
340  fullDomainMapVector.insert(fullDomainMapVector.end(), nodeDomainMapii.begin(), nodeDomainMapii.end());
341  if(bThyraDomainGIDs == false)
342  sort(fullDomainMapVector.begin(), fullDomainMapVector.end());
343 
344  curBlockId++; // increase block id index
345 
346  } // end SubFactoryManager environment
347 
348  // extract map index base from maps of blocked P
349  GO rangeIndexBase = originalTransferOp->getRangeMap()->getIndexBase();
350  GO domainIndexBase= originalTransferOp->getDomainMap()->getIndexBase();
351 
352  RCP<const Xpetra::MapExtractor<Scalar, LocalOrdinal, GlobalOrdinal, Node> > rangePMapExtractor = bOriginalTransferOp->getRangeMapExtractor(); // original map extractor
353  Teuchos::ArrayView<GO> fullRangeMapGIDs(fullRangeMapVector.size() ? &fullRangeMapVector[0] : 0,fullRangeMapVector.size());
354  Teuchos::RCP<const StridedMap> stridedRgFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(rangePMapExtractor->getFullMap());
355  Teuchos::RCP<const Map > fullRangeMap = Teuchos::null;
356  if(stridedRgFullMap != Teuchos::null) {
357  std::vector<size_t> stridedData = stridedRgFullMap->getStridingData();
358  fullRangeMap =
359  StridedMapFactory::Build(
360  originalTransferOp->getRangeMap()->lib(),
361  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
362  fullRangeMapGIDs,
363  rangeIndexBase,
364  stridedData,
365  originalTransferOp->getRangeMap()->getComm(),
366  stridedRgFullMap->getStridedBlockId(),
367  stridedRgFullMap->getOffset());
368  } else {
369  fullRangeMap =
370  MapFactory::Build(
371  originalTransferOp->getRangeMap()->lib(),
372  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
373  fullRangeMapGIDs,
374  rangeIndexBase,
375  originalTransferOp->getRangeMap()->getComm());
376  }
377 
378  RCP<const Xpetra::MapExtractor<Scalar, LocalOrdinal, GlobalOrdinal, Node> > domainAMapExtractor = bOriginalTransferOp->getDomainMapExtractor();
379  Teuchos::ArrayView<GO> fullDomainMapGIDs(fullDomainMapVector.size() ? &fullDomainMapVector[0] : 0,fullDomainMapVector.size());
380  Teuchos::RCP<const StridedMap> stridedDoFullMap = Teuchos::rcp_dynamic_cast<const StridedMap>(domainAMapExtractor->getFullMap());
381  Teuchos::RCP<const Map > fullDomainMap = Teuchos::null;
382  if(stridedDoFullMap != Teuchos::null) {
383  TEUCHOS_TEST_FOR_EXCEPTION(stridedDoFullMap==Teuchos::null, Exceptions::BadCast, "MueLu::BlockedPFactory::Build: full map in domain map extractor has no striding information! error.");
384  std::vector<size_t> stridedData2 = stridedDoFullMap->getStridingData();
385  fullDomainMap =
386  StridedMapFactory::Build(
387  originalTransferOp->getDomainMap()->lib(),
388  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
389  fullDomainMapGIDs,
390  domainIndexBase,
391  stridedData2,
392  originalTransferOp->getDomainMap()->getComm(),
393  stridedDoFullMap->getStridedBlockId(),
394  stridedDoFullMap->getOffset());
395  } else {
396 
397  fullDomainMap =
398  MapFactory::Build(
399  originalTransferOp->getDomainMap()->lib(),
400  Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),
401  fullDomainMapGIDs,
402  domainIndexBase,
403  originalTransferOp->getDomainMap()->getComm());
404  }
405 
406  // build map extractors
407  Teuchos::RCP<const MapExtractor> rebrangeMapExtractor =
408  MapExtractorFactory::Build(fullRangeMap, subBlockPRangeMaps, bThyraRangeGIDs);
409  Teuchos::RCP<const MapExtractor> rebdomainMapExtractor =
410  MapExtractorFactory::Build(fullDomainMap, subBlockPDomainMaps, bThyraDomainGIDs);
411 
412  Teuchos::RCP<BlockedCrsMatrix> bRebP = Teuchos::rcp(new BlockedCrsMatrix(rebrangeMapExtractor,rebdomainMapExtractor,10));
413  for(size_t i = 0; i<subBlockPRangeMaps.size(); i++) {
414  Teuchos::RCP<CrsMatrixWrap> crsOpii = Teuchos::rcp_dynamic_cast<CrsMatrixWrap>(subBlockRebP[i]);
415  TEUCHOS_TEST_FOR_EXCEPTION(crsOpii == Teuchos::null,Xpetra::Exceptions::BadCast, "MueLu::RebalanceBlockTransferFactory::Build: block P" << i << " is not of type CrsMatrixWrap.");
416  bRebP->setMatrix(i,i,crsOpii);
417  }
418  bRebP->fillComplete();
419 
420  Set(coarseLevel, "P", Teuchos::rcp_dynamic_cast<Matrix>(bRebP));
421 } // Build
422 
423 } // namespace MueLu
424 
425 #endif /* HAVE_MUELU_EXPERIMENTAL */
426 #endif /* MUELU_REBALANCEBLOCKINTERPOLATIONFACTORY_DEF_HPP_ */
Exception indicating invalid cast attempted.
#define MueLu_maxAll(rcpComm, in, out)
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Timer to be used in factories. Similar to Monitor but with additional timers.
Namespace for MueLu classes and methods.
Print statistics that do not involve significant additional computation.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
Exception throws to report errors in the internal logical of the program.
#define MUELU_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
An exception safe way to call the method &#39;Level::SetFactoryManager()&#39;.
void AddFactoryManager(RCP< const FactoryManagerBase > FactManager)
Add a factory manager.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()