46 #ifndef MUELU_UNCOUPLEDAGGREGATIONFACTORY_KOKKOS_DEF_HPP_
47 #define MUELU_UNCOUPLEDAGGREGATIONFACTORY_KOKKOS_DEF_HPP_
49 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
53 #include <Xpetra_Map.hpp>
55 #include <Xpetra_MultiVectorFactory.hpp>
60 #include "MueLu_OnePtAggregationAlgorithm_kokkos.hpp"
61 #include "MueLu_PreserveDirichletAggregationAlgorithm_kokkos.hpp"
62 #include "MueLu_IsolatedNodeAggregationAlgorithm_kokkos.hpp"
64 #include "MueLu_AggregationPhase1Algorithm_kokkos.hpp"
65 #include "MueLu_AggregationPhase2aAlgorithm_kokkos.hpp"
66 #include "MueLu_AggregationPhase2bAlgorithm_kokkos.hpp"
67 #include "MueLu_AggregationPhase3Algorithm_kokkos.hpp"
70 #include "MueLu_LWGraph_kokkos.hpp"
71 #include "MueLu_Aggregates_kokkos.hpp"
74 #include "MueLu_AmalgamationInfo.hpp"
75 #include "MueLu_Utilities.hpp"
77 #include "KokkosGraph_Distance2ColorHandle.hpp"
78 #include "KokkosGraph_Distance2Color.hpp"
82 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
83 UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::UncoupledAggregationFactory_kokkos()
84 : bDefinitionPhase_(true)
87 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
88 RCP<const ParameterList> UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::GetValidParameterList()
const {
89 RCP<ParameterList> validParamList =
rcp(
new ParameterList());
95 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
100 validParamList->getEntry(
"aggregation: ordering").setValidator(
101 rcp(
new validatorType(Teuchos::tuple<std::string>(
"natural",
"graph",
"random"),
"aggregation: ordering")));
110 #undef SET_VALID_ENTRY
113 validParamList->set< RCP<const FactoryBase> >(
"Graph", null,
"Generating factory of the graph");
114 validParamList->set< RCP<const FactoryBase> >(
"DofsPerNode", null,
"Generating factory for variable \'DofsPerNode\', usually the same as for \'Graph\'");
117 validParamList->set< std::string > (
"OnePt aggregate map name",
"",
"Name of input map for single node aggregates. (default='')");
118 validParamList->set< std::string > (
"OnePt aggregate map factory",
"",
"Generating factory of (DOF) map for single node aggregates.");
121 return validParamList;
124 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 void UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::DeclareInput(Level& currentLevel)
const {
126 Input(currentLevel,
"Graph");
127 Input(currentLevel,
"DofsPerNode");
129 const ParameterList& pL = GetParameterList();
132 std::string mapOnePtName = pL.get<std::string>(
"OnePt aggregate map name");
133 if (mapOnePtName.length() > 0) {
134 std::string mapOnePtFactName = pL.get<std::string>(
"OnePt aggregate map factory");
135 if (mapOnePtFactName ==
"" || mapOnePtFactName ==
"NoFactory") {
136 currentLevel.DeclareInput(mapOnePtName, NoFactory::get());
138 RCP<const FactoryBase> mapOnePtFact = GetFactory(mapOnePtFactName);
139 currentLevel.DeclareInput(mapOnePtName, mapOnePtFact.get());
144 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
145 void UncoupledAggregationFactory_kokkos<LocalOrdinal, GlobalOrdinal, Node>::Build(Level ¤tLevel)
const {
146 FactoryMonitor m(*
this,
"Build", currentLevel);
148 ParameterList pL = GetParameterList();
149 bDefinitionPhase_ =
false;
151 if (pL.get<
int>(
"aggregation: max agg size") == -1)
152 pL.set(
"aggregation: max agg size", INT_MAX);
155 RCP<const FactoryBase> graphFact = GetFactory(
"Graph");
159 algos_.push_back(
rcp(
new PreserveDirichletAggregationAlgorithm_kokkos(graphFact)));
160 if (pL.get<
bool>(
"aggregation: allow user-specified singletons") ==
true) algos_.push_back(
rcp(
new OnePtAggregationAlgorithm_kokkos (graphFact)));
161 if (pL.get<
bool>(
"aggregation: enable phase 1" ) ==
true) algos_.push_back(
rcp(
new AggregationPhase1Algorithm_kokkos (graphFact)));
162 if (pL.get<
bool>(
"aggregation: enable phase 2a") ==
true) algos_.push_back(
rcp(
new AggregationPhase2aAlgorithm_kokkos (graphFact)));
163 if (pL.get<
bool>(
"aggregation: enable phase 2b") ==
true) algos_.push_back(
rcp(
new AggregationPhase2bAlgorithm_kokkos (graphFact)));
164 if (pL.get<
bool>(
"aggregation: enable phase 3" ) ==
true) algos_.push_back(
rcp(
new AggregationPhase3Algorithm_kokkos (graphFact)));
166 std::string mapOnePtName = pL.get<std::string>(
"OnePt aggregate map name");
167 RCP<Map> OnePtMap = Teuchos::null;
168 if (mapOnePtName.length()) {
169 std::string mapOnePtFactName = pL.get<std::string>(
"OnePt aggregate map factory");
170 if (mapOnePtFactName ==
"" || mapOnePtFactName ==
"NoFactory") {
171 OnePtMap = currentLevel.Get<RCP<Map> >(mapOnePtName, NoFactory::get());
173 RCP<const FactoryBase> mapOnePtFact = GetFactory(mapOnePtFactName);
174 OnePtMap = currentLevel.Get<RCP<Map> >(mapOnePtName, mapOnePtFact.get());
178 RCP<const LWGraph_kokkos> graph = Get< RCP<LWGraph_kokkos> >(currentLevel,
"Graph");
181 RCP<Aggregates_kokkos> aggregates =
rcp(
new Aggregates_kokkos(*graph));
182 aggregates->setObjectLabel(
"UC");
184 const LO numRows = graph->GetNodeNumVertices();
208 LO nDofsPerNode = Get<LO>(currentLevel,
"DofsPerNode");
209 GO indexBase = graph->GetDomainMap()->getIndexBase();
210 if (OnePtMap != Teuchos::null) {
212 = Kokkos::create_mirror_view(aggStat);
215 for (
LO i = 0; i < numRows; i++) {
217 GO grid = (graph->GetDomainMap()->getGlobalElement(i)-indexBase) * nDofsPerNode + indexBase;
219 for (
LO kr = 0; kr < nDofsPerNode; kr++)
220 if (OnePtMap->isNodeGlobalElement(grid + kr))
221 aggStatHost(i) =
ONEPT;
228 const RCP<const Teuchos::Comm<int> > comm = graph->GetComm();
229 GO numGlobalRows = 0;
234 SubFactoryMonitor sfm(*
this,
"Algo \"Graph Coloring\"", currentLevel);
241 using graph_t =
typename LWGraph_kokkos::local_graph_type;
242 using KernelHandle = KokkosKernels::Experimental::
243 KokkosKernelsHandle<
typename graph_t::row_map_type::value_type,
244 typename graph_t::entries_type::value_type,
245 typename graph_t::entries_type::value_type,
246 typename graph_t::device_type::execution_space,
247 typename graph_t::device_type::memory_space,
248 typename graph_t::device_type::memory_space>;
251 kh.create_distance2_graph_coloring_handle();
254 auto coloringHandle = kh.get_distance2_graph_coloring_handle();
265 coloringHandle->set_algorithm( KokkosGraph::COLORING_D2_SERIAL );
268 typename graph_t::row_map_type aRowptrs = graph->getRowPtrs();
269 typename graph_t::entries_type aColinds = graph->getEntries();
273 KokkosGraph::Experimental::graph_compute_distance2_color(&kh, numRows, numRows,
278 aggregates->SetGraphColors(coloringHandle->get_vertex_colors());
279 aggregates->SetGraphNumColors(static_cast<LO>(coloringHandle->get_num_colors()));
282 kh.destroy_distance2_graph_coloring_handle();
285 LO numNonAggregatedNodes = numRows;
286 GO numGlobalAggregatedPrev = 0, numGlobalAggsPrev = 0;
287 for (
size_t a = 0; a < algos_.size(); a++) {
288 std::string phase = algos_[a]->description();
289 SubFactoryMonitor sfm(*
this,
"Algo \"" + phase +
"\"", currentLevel);
291 int oldRank = algos_[a]->SetProcRankVerbose(this->GetProcRankVerbose());
292 algos_[a]->BuildAggregates(pL, *graph, *aggregates, aggStat, numNonAggregatedNodes);
293 algos_[a]->SetProcRankVerbose(oldRank);
296 GO numLocalAggregated = numRows - numNonAggregatedNodes, numGlobalAggregated = 0;
297 GO numLocalAggs = aggregates->GetNumAggregates(), numGlobalAggs = 0;
298 MueLu_sumAll(comm, numLocalAggregated, numGlobalAggregated);
301 double aggPercent = 100*as<double>(numGlobalAggregated)/as<double>(numGlobalRows);
302 if (aggPercent > 99.99 && aggPercent < 100.00) {
309 GetOStream(
Statistics1) <<
" aggregated : " << (numGlobalAggregated - numGlobalAggregatedPrev) <<
" (phase), " << std::fixed
310 << std::setprecision(2) << numGlobalAggregated <<
"/" << numGlobalRows <<
" [" << aggPercent <<
"%] (total)\n"
311 <<
" remaining : " << numGlobalRows - numGlobalAggregated <<
"\n"
312 <<
" aggregates : " << numGlobalAggs-numGlobalAggsPrev <<
" (phase), " << numGlobalAggs <<
" (total)" << std::endl;
313 numGlobalAggregatedPrev = numGlobalAggregated;
314 numGlobalAggsPrev = numGlobalAggs;
318 TEUCHOS_TEST_FOR_EXCEPTION(numNonAggregatedNodes, Exceptions::RuntimeError,
"MueLu::UncoupledAggregationFactory::Build: Leftover nodes found! Error!");
320 aggregates->AggregatesCrossProcessors(
false);
321 aggregates->ComputeAggregateSizes(
true);
323 Set(currentLevel,
"Aggregates", aggregates);
325 GetOStream(
Statistics1) << aggregates->description() << std::endl;
330 #endif // HAVE_MUELU_KOKKOS_REFACTOR
#define MueLu_sumAll(rcpComm, in, out)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void deep_copy(const View< DT, DP...> &dst, typename ViewTraits< DT, DP...>::const_value_type &value, typename std::enable_if< std::is_same< typename ViewTraits< DT, DP...>::specialize, void >::value >::type *=0)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
#define SET_VALID_ENTRY(name)