MueLu Version of the Day
MueLu_RebalanceAcFactory_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_REBALANCEACFACTORY_DEF_HPP
47#define MUELU_REBALANCEACFACTORY_DEF_HPP
48
49#include <Xpetra_Matrix.hpp>
50#include <Xpetra_CrsMatrix.hpp>
51#include <Xpetra_CrsMatrixWrap.hpp>
52#include <Xpetra_MatrixFactory.hpp>
53
55
56#include "MueLu_MasterList.hpp"
57#include "MueLu_Monitor.hpp"
58#include "MueLu_PerfUtils.hpp"
59#include "MueLu_RAPFactory.hpp"
60
61namespace MueLu {
62
63 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
65 RCP<ParameterList> validParamList = rcp(new ParameterList());
66
67#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
68 SET_VALID_ENTRY("repartition: use subcommunicators");
69#undef SET_VALID_ENTRY
70
71 validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A for rebalancing");
72 validParamList->set< RCP<const FactoryBase> >("Importer", Teuchos::null, "Generating factory of the importer");
73
74 return validParamList;
75 }
76
77 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
79 Input(coarseLevel, "A");
80 Input(coarseLevel, "Importer");
81 }
82
83 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
85 FactoryMonitor m(*this, "Computing Ac", coarseLevel);
86
87 RCP<Matrix> originalAc = Get< RCP<Matrix> >(coarseLevel, "A");
88 RCP<const Import> rebalanceImporter = Get< RCP<const Import> >(coarseLevel, "Importer");
89
90 if (rebalanceImporter != Teuchos::null) {
91 RCP<Matrix> rebalancedAc;
92 {
93 SubFactoryMonitor subM(*this, "Rebalancing existing Ac", coarseLevel);
94 RCP<const Map> targetMap = rebalanceImporter->getTargetMap();
95
96 const ParameterList & pL = GetParameterList();
97
98 ParameterList XpetraList;
99 if (pL.get<bool>("repartition: use subcommunicators") == true) {
100 GetOStream(Runtime0) << "Replacing maps with a subcommunicator" << std::endl;
101 XpetraList.set("Restrict Communicator",true);
102 }
103 // NOTE: If the communicator is restricted away, Build returns Teuchos::null.
104 XpetraList.set("Timer Label","MueLu::RebalanceAc-" + Teuchos::toString(coarseLevel.GetLevelID()));
105 rebalancedAc = MatrixFactory::Build(originalAc, *rebalanceImporter, *rebalanceImporter, targetMap, targetMap, rcp(&XpetraList,false));
106
107 if (!rebalancedAc.is_null()) {
108 rebalancedAc->SetFixedBlockSize(originalAc->GetFixedBlockSize());
109 std::ostringstream oss; oss << "A_" << coarseLevel.GetLevelID(); rebalancedAc->setObjectLabel(oss.str());
110 }
111 Set(coarseLevel, "A", rebalancedAc);
112 }
113 if (!rebalancedAc.is_null() && IsPrint(Statistics2)) {
114 int oldRank = SetProcRankVerbose(rebalancedAc->getRowMap()->getComm()->getRank());
115
116 RCP<ParameterList> params = rcp(new ParameterList());
117 params->set("printLoadBalancingInfo", true);
118 params->set("printCommInfo", true);
119 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*rebalancedAc, "Ac (rebalanced)", params);
120
121 SetProcRankVerbose(oldRank);
122 }
123
124 } else {
125 // Ac already built by the load balancing process and no load balancing needed
126 GetOStream(Runtime1) << "No rebalancing" << std::endl;
127 Set(coarseLevel, "A", originalAc);
128 }
129
130 if (rebalanceFacts_.begin() != rebalanceFacts_.end()) {
131 SubFactoryMonitor m2(*this, "Rebalance additional data", coarseLevel);
132
133 // call Build of all user-given transfer factories
134 for (std::vector<RCP<const FactoryBase> >::const_iterator it = rebalanceFacts_.begin(); it != rebalanceFacts_.end(); ++it) {
135 GetOStream(Runtime0) << "RebalanceAc: call rebalance factory " << (*it).get() << ": " << (*it)->description() << std::endl;
136 (*it)->CallBuild(coarseLevel);
137 }
138 }
139 } //Build()
140
141
142 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
144
145 /*TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<const TwoLevelFactoryBase>(factory) == Teuchos::null, Exceptions::BadCast,
146 "MueLu::RAPFactory::AddTransferFactory: Transfer factory is not derived from TwoLevelFactoryBase. "
147 "This is very strange. (Note: you can remove this exception if there's a good reason for)");
148 TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_, Exceptions::RuntimeError, "MueLu::RAPFactory::AddTransferFactory: Factory is being added after we have already declared input");*/
149 rebalanceFacts_.push_back(factory);
150 } //AddRebalanceFactory()
151
152} //namespace MueLu
153
154#endif // MUELU_REBALANCEACFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
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
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void AddRebalanceFactory(const RCP< const FactoryBase > &factory)
Add rebalancing factory in the end of list of rebalancing factories in RebalanceAcFactory.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
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.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.