MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_UncoupledAggregationFactory_kokkos_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_UNCOUPLEDAGGREGATIONFACTORY_KOKKOS_DEF_HPP_
47 #define MUELU_UNCOUPLEDAGGREGATIONFACTORY_KOKKOS_DEF_HPP_
48 
49 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
50 
51 #include <climits>
52 
53 #include <Xpetra_Map.hpp>
54 #include <Xpetra_Vector.hpp>
56 #include <Xpetra_VectorFactory.hpp>
57 
59 
60 #include "MueLu_OnePtAggregationAlgorithm_kokkos.hpp"
61 #include "MueLu_PreserveDirichletAggregationAlgorithm_kokkos.hpp"
62 #include "MueLu_IsolatedNodeAggregationAlgorithm_kokkos.hpp"
63 
64 #include "MueLu_AggregationPhase1Algorithm_kokkos.hpp"
65 #include "MueLu_AggregationPhase2aAlgorithm_kokkos.hpp"
66 #include "MueLu_AggregationPhase2bAlgorithm_kokkos.hpp"
67 #include "MueLu_AggregationPhase3Algorithm_kokkos.hpp"
68 
69 #include "MueLu_Level.hpp"
70 #include "MueLu_LWGraph_kokkos.hpp"
71 #include "MueLu_Aggregates_kokkos.hpp"
72 #include "MueLu_MasterList.hpp"
73 #include "MueLu_Monitor.hpp"
74 #include "MueLu_AmalgamationInfo.hpp"
75 #include "MueLu_Utilities.hpp" // for sum_all and similar stuff...
76 
77 namespace MueLu {
78 
79  template <class LocalOrdinal, class GlobalOrdinal, class Node>
80  UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::UncoupledAggregationFactory_kokkos()
81  : bDefinitionPhase_(true)
82  { }
83 
84  template <class LocalOrdinal, class GlobalOrdinal, class Node>
85  RCP<const ParameterList> UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::GetValidParameterList() const {
86  RCP<ParameterList> validParamList = rcp(new ParameterList());
87 
88  // Aggregation parameters (used in aggregation algorithms)
89  // TODO introduce local member function for each aggregation algorithm such that each aggregation algorithm can define its own parameters
90 
92 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
93  SET_VALID_ENTRY("aggregation: max agg size");
94  SET_VALID_ENTRY("aggregation: min agg size");
95  SET_VALID_ENTRY("aggregation: max selected neighbors");
96  SET_VALID_ENTRY("aggregation: ordering");
97  validParamList->getEntry("aggregation: ordering").setValidator(
98  rcp(new validatorType(Teuchos::tuple<std::string>("natural", "graph", "random"), "aggregation: ordering")));
99  SET_VALID_ENTRY("aggregation: enable phase 1");
100  SET_VALID_ENTRY("aggregation: phase 1 algorithm");
101  SET_VALID_ENTRY("aggregation: enable phase 2a");
102  SET_VALID_ENTRY("aggregation: enable phase 2b");
103  SET_VALID_ENTRY("aggregation: enable phase 3");
104  SET_VALID_ENTRY("aggregation: preserve Dirichlet points");
105  SET_VALID_ENTRY("aggregation: allow user-specified singletons");
106 #undef SET_VALID_ENTRY
107 
108  // general variables needed in AggregationFactory
109  validParamList->set< RCP<const FactoryBase> >("Graph", null, "Generating factory of the graph");
110  validParamList->set< RCP<const FactoryBase> >("DofsPerNode", null, "Generating factory for variable \'DofsPerNode\', usually the same as for \'Graph\'");
111 
112  // special variables necessary for OnePtAggregationAlgorithm
113  validParamList->set< std::string > ("OnePt aggregate map name", "", "Name of input map for single node aggregates. (default='')");
114  validParamList->set< std::string > ("OnePt aggregate map factory", "", "Generating factory of (DOF) map for single node aggregates.");
115  //validParamList->set< RCP<const FactoryBase> >("OnePt aggregate map factory", NoFactory::getRCP(), "Generating factory of (DOF) map for single node aggregates.");
116 
117  return validParamList;
118  }
119 
120  template <class LocalOrdinal, class GlobalOrdinal, class Node>
121  void UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::DeclareInput(Level& currentLevel) const {
122  Input(currentLevel, "Graph");
123  Input(currentLevel, "DofsPerNode");
124 
125  const ParameterList& pL = GetParameterList();
126 
127  // request special data necessary for OnePtAggregationAlgorithm
128  std::string mapOnePtName = pL.get<std::string>("OnePt aggregate map name");
129  if (mapOnePtName.length() > 0) {
130  std::string mapOnePtFactName = pL.get<std::string>("OnePt aggregate map factory");
131  if (mapOnePtFactName == "" || mapOnePtFactName == "NoFactory") {
132  currentLevel.DeclareInput(mapOnePtName, NoFactory::get());
133  } else {
134  RCP<const FactoryBase> mapOnePtFact = GetFactory(mapOnePtFactName);
135  currentLevel.DeclareInput(mapOnePtName, mapOnePtFact.get());
136  }
137  }
138  }
139 
140  template <class LocalOrdinal, class GlobalOrdinal, class Node>
141  void UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::Build(Level &currentLevel) const {
142  FactoryMonitor m(*this, "Build", currentLevel);
143 
144  ParameterList pL = GetParameterList();
145  bDefinitionPhase_ = false; // definition phase is finished, now all aggregation algorithm information is fixed
146 
147  if (pL.get<int>("aggregation: max agg size") == -1)
148  pL.set("aggregation: max agg size", INT_MAX);
149 
150  // define aggregation algorithms
151  RCP<const FactoryBase> graphFact = GetFactory("Graph");
152 
153  // TODO Can we keep different aggregation algorithms over more Build calls?
154  algos_.clear();
155  algos_.push_back(rcp(new PreserveDirichletAggregationAlgorithm_kokkos(graphFact)));
156  if (pL.get<bool>("aggregation: allow user-specified singletons") == true) algos_.push_back(rcp(new OnePtAggregationAlgorithm_kokkos (graphFact)));
157  if (pL.get<bool>("aggregation: enable phase 1" ) == true) algos_.push_back(rcp(new AggregationPhase1Algorithm_kokkos (graphFact)));
158  if (pL.get<bool>("aggregation: enable phase 2a") == true) algos_.push_back(rcp(new AggregationPhase2aAlgorithm_kokkos (graphFact)));
159  if (pL.get<bool>("aggregation: enable phase 2b") == true) algos_.push_back(rcp(new AggregationPhase2bAlgorithm_kokkos (graphFact)));
160  if (pL.get<bool>("aggregation: enable phase 3" ) == true) algos_.push_back(rcp(new AggregationPhase3Algorithm_kokkos (graphFact)));
161 
162  std::string mapOnePtName = pL.get<std::string>("OnePt aggregate map name");
163  RCP<Map> OnePtMap = Teuchos::null;
164  if (mapOnePtName.length()) {
165  std::string mapOnePtFactName = pL.get<std::string>("OnePt aggregate map factory");
166  if (mapOnePtFactName == "" || mapOnePtFactName == "NoFactory") {
167  OnePtMap = currentLevel.Get<RCP<Map> >(mapOnePtName, NoFactory::get());
168  } else {
169  RCP<const FactoryBase> mapOnePtFact = GetFactory(mapOnePtFactName);
170  OnePtMap = currentLevel.Get<RCP<Map> >(mapOnePtName, mapOnePtFact.get());
171  }
172  }
173 
174  RCP<const LWGraph_kokkos> graph = Get< RCP<LWGraph_kokkos> >(currentLevel, "Graph");
175 
176  // Build
177  RCP<Aggregates_kokkos> aggregates = rcp(new Aggregates_kokkos(*graph));
178  aggregates->setObjectLabel("UC");
179 
180  const LO numRows = graph->GetNodeNumVertices();
181 
182  // construct aggStat information
183  std::vector<unsigned> aggStat(numRows, READY);
184 
185  // TODO
186  //ArrayRCP<const bool> dirichletBoundaryMap = graph->GetBoundaryNodeMap();
187  ArrayRCP<const bool> dirichletBoundaryMap;
188 
189  if (dirichletBoundaryMap != Teuchos::null)
190  for (LO i = 0; i < numRows; i++)
191  if (dirichletBoundaryMap[i] == true)
192  aggStat[i] = BOUNDARY;
193 
194  LO nDofsPerNode = Get<LO>(currentLevel, "DofsPerNode");
195  GO indexBase = graph->GetDomainMap()->getIndexBase();
196  if (OnePtMap != Teuchos::null) {
197  for (LO i = 0; i < numRows; i++) {
198  // reconstruct global row id (FIXME only works for contiguous maps)
199  GO grid = (graph->GetDomainMap()->getGlobalElement(i)-indexBase) * nDofsPerNode + indexBase;
200 
201  for (LO kr = 0; kr < nDofsPerNode; kr++)
202  if (OnePtMap->isNodeGlobalElement(grid + kr))
203  aggStat[i] = ONEPT;
204  }
205  }
206 
207 
208  const RCP<const Teuchos::Comm<int> > comm = graph->GetComm();
209  GO numGlobalRows = 0;
210  if (IsPrint(Statistics1))
211  MueLu_sumAll(comm, as<GO>(numRows), numGlobalRows);
212 
213  LO numNonAggregatedNodes = numRows;
214  GO numGlobalAggregatedPrev = 0, numGlobalAggsPrev = 0;
215  for (size_t a = 0; a < algos_.size(); a++) {
216  std::string phase = algos_[a]->description();
217  SubFactoryMonitor sfm(*this, "Algo \"" + phase + "\"", currentLevel);
218 
219  int oldRank = algos_[a]->SetProcRankVerbose(this->GetProcRankVerbose());
220  algos_[a]->BuildAggregates(pL, *graph, *aggregates, aggStat, numNonAggregatedNodes);
221  algos_[a]->SetProcRankVerbose(oldRank);
222 
223  if (IsPrint(Statistics1)) {
224  GO numLocalAggregated = numRows - numNonAggregatedNodes, numGlobalAggregated = 0;
225  GO numLocalAggs = aggregates->GetNumAggregates(), numGlobalAggs = 0;
226  MueLu_sumAll(comm, numLocalAggregated, numGlobalAggregated);
227  MueLu_sumAll(comm, numLocalAggs, numGlobalAggs);
228 
229  double aggPercent = 100*as<double>(numGlobalAggregated)/as<double>(numGlobalRows);
230  if (aggPercent > 99.99 && aggPercent < 100.00) {
231  // Due to round off (for instance, for 140465733/140466897), we could
232  // get 100.00% display even if there are some remaining nodes. This
233  // is bad from the users point of view. It is much better to change
234  // it to display 99.99%.
235  aggPercent = 99.99;
236  }
237  GetOStream(Statistics1) << " aggregated : " << (numGlobalAggregated - numGlobalAggregatedPrev) << " (phase), " << std::fixed
238  << std::setprecision(2) << numGlobalAggregated << "/" << numGlobalRows << " [" << aggPercent << "%] (total)\n"
239  << " remaining : " << numGlobalRows - numGlobalAggregated << "\n"
240  << " aggregates : " << numGlobalAggs-numGlobalAggsPrev << " (phase), " << numGlobalAggs << " (total)" << std::endl;
241  numGlobalAggregatedPrev = numGlobalAggregated;
242  numGlobalAggsPrev = numGlobalAggs;
243  }
244  }
245 
246  TEUCHOS_TEST_FOR_EXCEPTION(numNonAggregatedNodes, Exceptions::RuntimeError, "MueLu::UncoupledAggregationFactory::Build: Leftover nodes found! Error!");
247 
248  aggregates->AggregatesCrossProcessors(false);
249  aggregates->ComputeAggregateSizes(true/*forceRecompute*/);
250 
251  Set(currentLevel, "Aggregates", aggregates);
252 
253  GetOStream(Statistics1) << aggregates->description() << std::endl;
254  }
255 
256 } //namespace MueLu
257 
258 #endif // HAVE_MUELU_KOKKOS_REFACTOR
259 #endif /* MUELU_UNCOUPLEDAGGREGATIONFACTORY_DEF_HPP_ */
#define MueLu_sumAll(rcpComm, in, out)
GlobalOrdinal GO
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Print more statistics.
LocalOrdinal LO
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
#define SET_VALID_ENTRY(name)