10 #ifndef MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
11 #define MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
14 #include "Xpetra_MultiVectorFactory.hpp"
15 #include "Xpetra_MapFactory.hpp"
18 #include "MueLu_Aggregates.hpp"
19 #include "MueLu_AmalgamationFactory.hpp"
21 #include "MueLu_Utilities.hpp"
28 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
31 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
34 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
41 validParamList->
set<
bool>(
"structured aggregation",
false,
"Flag specifying that the geometric data is transferred for StructuredAggregationFactory");
42 validParamList->
set<
bool>(
"aggregation coupled",
false,
"Flag specifying if the aggregation algorithm was used in coupled mode.");
43 validParamList->
set<
bool>(
"Geometric",
false,
"Flag specifying that the coordinates are transferred for GeneralGeometricPFactory");
44 validParamList->
set<
RCP<const FactoryBase> >(
"coarseCoordinates", Teuchos::null,
"Factory for coarse coordinates generation");
45 validParamList->
set<
RCP<const FactoryBase> >(
"gCoarseNodesPerDim", Teuchos::null,
"Factory providing the global number of nodes per spatial dimensions of the mesh");
46 validParamList->
set<
RCP<const FactoryBase> >(
"lCoarseNodesPerDim", Teuchos::null,
"Factory providing the local number of nodes per spatial dimensions of the mesh");
47 validParamList->
set<
RCP<const FactoryBase> >(
"numDimensions", Teuchos::null,
"Factory providing the number of spatial dimensions of the mesh");
48 validParamList->
set<
int>(
"write start", -1,
"first level at which coordinates should be written to file");
49 validParamList->
set<
int>(
"write end", -1,
"last level at which coordinates should be written to file");
50 validParamList->
set<
bool>(
"hybrid aggregation",
false,
"Flag specifying that hybrid aggregation data is transfered for HybridAggregationFactory");
51 validParamList->
set<
RCP<const FactoryBase> >(
"aggregationRegionTypeCoarse", Teuchos::null,
"Factory indicating what aggregation type is to be used on the coarse level of the region");
52 validParamList->
set<
bool>(
"interface aggregation",
false,
"Flag specifying that interface aggregation data is transfered for HybridAggregationFactory");
53 validParamList->
set<
RCP<const FactoryBase> >(
"coarseInterfacesDimensions", Teuchos::null,
"Factory providing coarseInterfacesDimensions");
54 validParamList->
set<
RCP<const FactoryBase> >(
"nodeOnCoarseInterface", Teuchos::null,
"Factory providing nodeOnCoarseInterface");
56 return validParamList;
59 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
61 static bool isAvailableCoords =
false;
64 if (pL.
get<
bool>(
"structured aggregation") ==
true) {
65 if (pL.
get<
bool>(
"aggregation coupled") ==
true) {
66 Input(fineLevel,
"gCoarseNodesPerDim");
68 Input(fineLevel,
"lCoarseNodesPerDim");
69 Input(fineLevel,
"numDimensions");
70 }
else if (pL.
get<
bool>(
"Geometric") ==
true) {
71 Input(coarseLevel,
"coarseCoordinates");
72 Input(coarseLevel,
"gCoarseNodesPerDim");
73 Input(coarseLevel,
"lCoarseNodesPerDim");
74 }
else if (pL.
get<
bool>(
"hybrid aggregation") ==
true) {
75 Input(fineLevel,
"aggregationRegionTypeCoarse");
76 Input(fineLevel,
"lCoarseNodesPerDim");
77 Input(fineLevel,
"numDimensions");
78 if (pL.
get<
bool>(
"interface aggregation") ==
true) {
79 Input(fineLevel,
"coarseInterfacesDimensions");
80 Input(fineLevel,
"nodeOnCoarseInterface");
84 isAvailableCoords = coarseLevel.
IsAvailable(
"Coordinates",
this);
86 if (isAvailableCoords ==
false) {
87 Input(fineLevel,
"Coordinates");
88 Input(fineLevel,
"Aggregates");
89 Input(fineLevel,
"CoarseMap");
94 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
100 GetOStream(
Runtime0) <<
"Transferring coordinates" << std::endl;
110 if (pL.
get<
bool>(
"hybrid aggregation") ==
true) {
111 std::string regionType = Get<std::string>(fineLevel,
"aggregationRegionTypeCoarse");
112 numDimensions = Get<int>(fineLevel,
"numDimensions");
113 lCoarseNodesPerDir = Get<Array<LO> >(fineLevel,
"lCoarseNodesPerDim");
114 Set<std::string>(coarseLevel,
"aggregationRegionType", regionType);
115 Set<int>(coarseLevel,
"numDimensions", numDimensions);
116 Set<Array<LO> >(coarseLevel,
"lNodesPerDim", lCoarseNodesPerDir);
118 if ((pL.
get<
bool>(
"interface aggregation") ==
true) && (regionType ==
"uncoupled")) {
119 Array<LO> coarseInterfacesDimensions = Get<Array<LO> >(fineLevel,
"coarseInterfacesDimensions");
120 Array<LO> nodeOnCoarseInterface = Get<Array<LO> >(fineLevel,
"nodeOnCoarseInterface");
121 Set<Array<LO> >(coarseLevel,
"interfacesDimensions", coarseInterfacesDimensions);
122 Set<Array<LO> >(coarseLevel,
"nodeOnInterface", nodeOnCoarseInterface);
125 }
else if (pL.
get<
bool>(
"structured aggregation") ==
true) {
126 if (pL.
get<
bool>(
"aggregation coupled") ==
true) {
127 gCoarseNodesPerDir = Get<Array<GO> >(fineLevel,
"gCoarseNodesPerDim");
128 Set<Array<GO> >(coarseLevel,
"gNodesPerDim", gCoarseNodesPerDir);
130 lCoarseNodesPerDir = Get<Array<LO> >(fineLevel,
"lCoarseNodesPerDim");
131 Set<Array<LO> >(coarseLevel,
"lNodesPerDim", lCoarseNodesPerDir);
132 numDimensions = Get<int>(fineLevel,
"numDimensions");
133 Set<int>(coarseLevel,
"numDimensions", numDimensions);
135 }
else if (pL.
get<
bool>(
"Geometric") ==
true) {
136 coarseCoords = Get<RCP<xdMV> >(coarseLevel,
"coarseCoordinates");
137 gCoarseNodesPerDir = Get<Array<GO> >(coarseLevel,
"gCoarseNodesPerDim");
138 lCoarseNodesPerDir = Get<Array<LO> >(coarseLevel,
"lCoarseNodesPerDim");
139 Set<Array<GO> >(coarseLevel,
"gNodesPerDim", gCoarseNodesPerDir);
140 Set<Array<LO> >(coarseLevel,
"lNodesPerDim", lCoarseNodesPerDir);
142 Set<RCP<xdMV> >(coarseLevel,
"Coordinates", coarseCoords);
145 if (coarseLevel.
IsAvailable(
"Coordinates",
this)) {
146 GetOStream(
Runtime0) <<
"Reusing coordinates" << std::endl;
150 fineCoords = Get<RCP<xdMV> >(fineLevel,
"Coordinates");
151 RCP<const Map> coarseMap = Get<RCP<const Map> >(fineLevel,
"CoarseMap");
160 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
161 blkSize = rcp_dynamic_cast<
const StridedMap>(coarseMap)->getFixedBlockSize();
163 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
164 blkSize = rcp_dynamic_cast<
const StridedMap>(coarseMap)->getFixedBlockSize();
169 coarseCoordMap = coarseMap;
179 bool aggregatesCrossProcessors;
180 aggregates = Get<RCP<Aggregates> >(fineLevel,
"Aggregates");
181 aggregatesCrossProcessors = aggregates->AggregatesCrossProcessors();
185 if (aggregatesCrossProcessors) {
195 auto aggGraph = aggregates->GetGraph();
196 auto numAggs = aggGraph.numRows();
198 auto fineCoordsView = ghostedCoords->getDeviceLocalView(Xpetra::Access::ReadOnly);
199 auto coarseCoordsView = coarseCoords->getDeviceLocalView(Xpetra::Access::OverwriteAll);
205 const auto dim = ghostedCoords->getNumVectors();
208 for (
size_t j = 0; j < dim; j++) {
209 Kokkos::parallel_for(
210 "MueLu:CoordinatesTransferF:Build:coord", Kokkos::RangePolicy<local_ordinal_type, execution_space>(0, numAggs),
211 KOKKOS_LAMBDA(
const LO i) {
215 auto aggregate = aggGraph.rowConst(i);
218 for (
size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
219 sum += fineCoordsRandomView(aggregate(colID), j);
221 coarseCoordsView(i, j) = sum / aggregate.length;
226 Set<RCP<xdMV> >(coarseLevel,
"Coordinates", coarseCoords);
229 int writeStart = pL.
get<
int>(
"write start"), writeEnd = pL.
get<
int>(
"write end");
230 if (writeStart == 0 && fineLevel.
GetLevelID() == 0 && writeStart <= writeEnd) {
231 std::ostringstream buf;
233 std::string fileName =
"coordinates_before_rebalance_level_" + buf.str() +
".m";
237 std::ostringstream buf;
239 std::string fileName =
"coordinates_before_rebalance_level_" + buf.str() +
".m";
246 #endif // MUELU_COORDINATESTRANSFER_FACTORY_DEF_HPP
CoordinatesTransferFactory()
Constructor.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
One-liner description of what is happening.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual ~CoordinatesTransferFactory()
Destructor.
RequestMode GetRequestMode() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static void AmalgamateMap(const Map &sourceMap, const Matrix &A, RCP< const Map > &amalgamatedMap, Array< LO > &translation)
Method to create merged map for systems of PDEs.
int GetLevelID() const
Return level number.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.