MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_StructuredAggregationFactory_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_STRUCTUREDAGGREGATIONFACTORY_KOKKOS_DEF_HPP
47#define MUELU_STRUCTUREDAGGREGATIONFACTORY_KOKKOS_DEF_HPP
48
49// Xpetra includes
50#include <Xpetra_Map.hpp>
51#include <Xpetra_CrsGraph.hpp>
52
53// MueLu generic includes
54#include "MueLu_Level.hpp"
55#include "MueLu_MasterList.hpp"
56#include "MueLu_Monitor.hpp"
57
58// MueLu specific includes (kokkos version)
59#include "MueLu_LWGraph_kokkos.hpp"
60#include "MueLu_Aggregates.hpp"
61#include "MueLu_IndexManager_kokkos.hpp"
62#include "MueLu_AggregationStructuredAlgorithm_kokkos.hpp"
63
65
66namespace MueLu {
67
68 template <class LocalOrdinal, class GlobalOrdinal, class Node>
71
72 template <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("aggregation: preserve Dirichlet points");
79 SET_VALID_ENTRY("aggregation: allow user-specified singletons");
80 SET_VALID_ENTRY("aggregation: error on nodes with no on-rank neighbors");
81 SET_VALID_ENTRY("aggregation: phase3 avoid singletons");
82#undef SET_VALID_ENTRY
83
84 // general variables needed in StructuredAggregationFactory
85 validParamList->set<std::string> ("aggregation: output type", "Aggregates",
86 "Type of object holding the aggregation data: Aggregtes or CrsGraph");
87 validParamList->set<std::string> ("aggregation: coarsening rate", "{3}",
88 "Coarsening rate per spatial dimensions");
89 validParamList->set<int> ("aggregation: coarsening order", 0,
90 "The interpolation order used to construct grid transfer operators based off these aggregates.");
91 validParamList->set<RCP<const FactoryBase> >("Graph", Teuchos::null,
92 "Graph of the matrix after amalgamation but without dropping.");
93 validParamList->set<RCP<const FactoryBase> >("DofsPerNode", Teuchos::null,
94 "Number of degrees of freedom per mesh node, provided by the coalsce drop factory.");
95 validParamList->set<RCP<const FactoryBase> >("numDimensions", Teuchos::null,
96 "Number of spatial dimension provided by CoordinatesTransferFactory.");
97 validParamList->set<RCP<const FactoryBase> >("lNodesPerDim", Teuchos::null,
98 "Number of nodes per spatial dimmension provided by CoordinatesTransferFactory.");
99
100 return validParamList;
101 } // GetValidParameterList()
102
103 template <class LocalOrdinal, class GlobalOrdinal, class Node>
105 DeclareInput(Level& currentLevel) const {
106 Input(currentLevel, "Graph");
107 Input(currentLevel, "DofsPerNode");
108
109 // Request the local number of nodes per dimensions
110 if(currentLevel.GetLevelID() == 0) {
111 if(currentLevel.IsAvailable("numDimensions", NoFactory::get())) {
112 currentLevel.DeclareInput("numDimensions", NoFactory::get(), this);
113 } else {
114 TEUCHOS_TEST_FOR_EXCEPTION(currentLevel.IsAvailable("numDimensions", NoFactory::get()),
116 "numDimensions was not provided by the user on level0!");
117 }
118 if(currentLevel.IsAvailable("lNodesPerDim", NoFactory::get())) {
119 currentLevel.DeclareInput("lNodesPerDim", NoFactory::get(), this);
120 } else {
121 TEUCHOS_TEST_FOR_EXCEPTION(currentLevel.IsAvailable("lNodesPerDim", NoFactory::get()),
123 "lNodesPerDim was not provided by the user on level0!");
124 }
125 } else {
126 Input(currentLevel, "lNodesPerDim");
127 Input(currentLevel, "numDimensions");
128 }
129 } // DeclareInput()
130
131 template <class LocalOrdinal, class GlobalOrdinal, class Node>
133 Build(Level &currentLevel) const {
134 FactoryMonitor m(*this, "Build", currentLevel);
135
136 RCP<Teuchos::FancyOStream> out;
137 if(const char* dbg = std::getenv("MUELU_STRUCTUREDAGGREGATION_DEBUG")) {
138 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
139 out->setShowAllFrontMatter(false).setShowProcRank(true);
140 } else {
141 out = Teuchos::getFancyOStream(rcp(new Teuchos::oblackholestream()));
142 }
143
144 using device_type = typename LWGraph_kokkos::local_graph_type::device_type;
145 using execution_space = typename LWGraph_kokkos::local_graph_type::device_type::execution_space;
146 using memory_space = typename LWGraph_kokkos::local_graph_type::device_type::memory_space;
147
148 *out << "Entering structured aggregation" << std::endl;
149
150 ParameterList pL = GetParameterList();
151 bDefinitionPhase_ = false; // definition phase is finished, now all aggregation algorithm information is fixed
152
153 // General problem informations are gathered from data stored in the problem matix.
154 RCP<const LWGraph_kokkos> graph = Get<RCP<LWGraph_kokkos> >(currentLevel, "Graph");
155 RCP<const Map> fineMap = graph->GetDomainMap();
156 const int myRank = fineMap->getComm()->getRank();
157 const LO dofsPerNode = Get<LO>(currentLevel, "DofsPerNode");
158
159 // Since we want to operate on nodes and not dof, we need to modify the rowMap in order to
160 // obtain a nodeMap.
161 const int interpolationOrder = pL.get<int>("aggregation: coarsening order");
162 std::string outputType = pL.get<std::string>("aggregation: output type");
163 const bool outputAggregates = (outputType == "Aggregates" ? true : false);
164 Array<LO> lFineNodesPerDir(3);
165 int numDimensions;
166 if(currentLevel.GetLevelID() == 0) {
167 // On level 0, data is provided by applications and has no associated factory.
168 lFineNodesPerDir = currentLevel.Get<Array<LO> >("lNodesPerDim", NoFactory::get());
169 numDimensions = currentLevel.Get<int>("numDimensions", NoFactory::get());
170 } else {
171 // On level > 0, data is provided directly by generating factories.
172 lFineNodesPerDir = Get<Array<LO> >(currentLevel, "lNodesPerDim");
173 numDimensions = Get<int>(currentLevel, "numDimensions");
174 }
175
176
177 // First make sure that input parameters are set logically based on dimension
178 for(int dim = 0; dim < 3; ++dim) {
179 if(dim >= numDimensions) {
180 lFineNodesPerDir[dim] = 1;
181 }
182 }
183
184 // Get the coarsening rate
185 std::string coarseningRate = pL.get<std::string>("aggregation: coarsening rate");
186 Teuchos::Array<LO> coarseRate;
187 try {
188 coarseRate = Teuchos::fromStringToArray<LO>(coarseningRate);
189 } catch(const Teuchos::InvalidArrayStringRepresentation& e) {
190 GetOStream(Errors,-1) << " *** \"aggregation: coarsening rate\" must be a string convertible into an array! *** "
191 << std::endl;
192 throw e;
193 }
194 TEUCHOS_TEST_FOR_EXCEPTION((coarseRate.size() > 1) && (coarseRate.size() < numDimensions),
196 "\"aggregation: coarsening rate\" must have at least as many"
197 " components as the number of spatial dimensions in the problem.");
198
199 // Now that we have extracted info from the level, create the IndexManager
200 RCP<IndexManager_kokkos> geoData = rcp(new IndexManager_kokkos(numDimensions,
201 interpolationOrder, myRank,
202 lFineNodesPerDir,
203 coarseRate));
204
205 *out << "The index manager has now been built" << std::endl;
206 TEUCHOS_TEST_FOR_EXCEPTION(fineMap->getLocalNumElements()
207 != static_cast<size_t>(geoData->getNumLocalFineNodes()),
209 "The local number of elements in the graph's map is not equal to "
210 "the number of nodes given by: lNodesPerDim!");
211
212 // Now we are ready for the big loop over the fine node that will assign each
213 // node on the fine grid to an aggregate and a processor.
214 RCP<AggregationStructuredAlgorithm_kokkos> myStructuredAlgorithm
216
217 if(interpolationOrder == 0 && outputAggregates){
218 RCP<Aggregates> aggregates = rcp(new Aggregates(graph->GetDomainMap()));
219 aggregates->setObjectLabel("ST");
220 aggregates->SetIndexManagerKokkos(geoData);
221 aggregates->AggregatesCrossProcessors(false);
222 aggregates->SetNumAggregates(geoData->getNumCoarseNodes());
223
224 LO numNonAggregatedNodes = geoData->getNumLocalFineNodes();
225 Kokkos::View<unsigned*, device_type> aggStat("aggStat", numNonAggregatedNodes);
226 Kokkos::parallel_for("StructuredAggregation: initialize aggStat",
227 Kokkos::RangePolicy<execution_space>(0, numNonAggregatedNodes),
228 KOKKOS_LAMBDA(const LO nodeIdx) {aggStat(nodeIdx) = READY;});
229
230 myStructuredAlgorithm->BuildAggregates(pL, *graph, *aggregates, aggStat,
231 numNonAggregatedNodes);
232
233 *out << "numNonAggregatedNodes: " << numNonAggregatedNodes << std::endl;
234
235 TEUCHOS_TEST_FOR_EXCEPTION(numNonAggregatedNodes, Exceptions::RuntimeError,
236 "MueLu::StructuredAggregationFactory::Build: Leftover nodes found! Error!");
237 aggregates->ComputeAggregateSizes(true/*forceRecompute*/);
238 GetOStream(Statistics1) << aggregates->description() << std::endl;
239 Set(currentLevel, "Aggregates", aggregates);
240
241 } else {
242 // Create Coarse Data
243 RCP<CrsGraph> myGraph;
244 myStructuredAlgorithm->BuildGraph(*graph, geoData, dofsPerNode, myGraph);
245 Set(currentLevel, "prolongatorGraph", myGraph);
246 }
247
248 Set(currentLevel, "lCoarseNodesPerDim", geoData->getCoarseNodesPerDirArray());
249 Set(currentLevel, "indexManager", geoData);
250 Set(currentLevel, "structuredInterpolationOrder", interpolationOrder);
251 Set(currentLevel, "numDimensions", numDimensions);
252
253 } // Build()
254
255} //namespace MueLu
256
257#endif /* MUELU_STRUCTUREDAGGREGATIONFACTORY_KOKKOS_DEF_HPP */
#define SET_VALID_ENTRY(name)
Container class for aggregation information.
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.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
Container class for mesh layout and indices calculation.
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput().
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)....
static const NoFactory * get()
virtual const Teuchos::ParameterList & GetParameterList() const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
Namespace for MueLu classes and methods.
@ Statistics1
Print more statistics.