MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_MultiVectorTransferFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
11#define MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
12
14#include "Tpetra_Access.hpp"
15#include "Xpetra_MultiVectorFactory.hpp"
16
17#include "MueLu_Aggregates.hpp"
18#include "MueLu_AmalgamationInfo.hpp"
19#include "MueLu_AmalgamationFactory.hpp"
20#include "MueLu_Level.hpp"
21#include "MueLu_UncoupledAggregationFactory.hpp"
22#include "MueLu_Monitor.hpp"
23
24namespace MueLu {
25
26template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28 RCP<ParameterList> validParamList = rcp(new ParameterList());
29
30 validParamList->set<std::string>("Vector name", "undefined", "Name of the vector that will be transferred on the coarse grid (level key)"); // TODO: how to set a validator without default value?
31 validParamList->set<std::string>("Transfer name", "R", "Name of the operator that will be used to transfer data");
32 validParamList->set<bool>("Normalize", false, "If a row sum normalization should be applied to preserve the mean value of the vector.");
33 validParamList->set<RCP<const FactoryBase>>("Vector factory", Teuchos::null, "Factory of the vector");
34 validParamList->set<RCP<const FactoryBase>>("Transfer factory", Teuchos::null, "Factory of the transfer operator");
35 validParamList->set<RCP<const FactoryBase>>("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
36
37 return validParamList;
38}
39
40template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42 const ParameterList &pL = GetParameterList();
43 std::string vectorName = pL.get<std::string>("Vector name");
44 std::string transferName = pL.get<std::string>("Transfer name");
45
46 fineLevel.DeclareInput(vectorName, GetFactory("Vector factory").get(), this);
47 auto transferFact = GetFactory("Transfer factory");
48 const bool isUncoupledAggFact = !Teuchos::rcp_dynamic_cast<const UncoupledAggregationFactory>(transferFact).is_null();
49 if (isUncoupledAggFact) {
50 fineLevel.DeclareInput(transferName, transferFact.get(), this);
51 Input(fineLevel, "CoarseMap");
52 } else
53 coarseLevel.DeclareInput(transferName, transferFact.get(), this);
54}
55
56template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
58 FactoryMonitor m(*this, "Build", coarseLevel);
59
60 const ParameterList &pL = GetParameterList();
61 const std::string vectorName = pL.get<std::string>("Vector name");
62 const std::string transferName = pL.get<std::string>("Transfer name");
63 const bool normalize = pL.get<bool>("Normalize");
64
65 auto transferFact = GetFactory("Transfer factory");
66 const bool isUncoupledAggFact = !Teuchos::rcp_dynamic_cast<const UncoupledAggregationFactory>(transferFact).is_null();
67
68 GetOStream(Runtime0) << "Transferring multivector \"" << vectorName << "\" using operator \"" << transferName << "\"" << std::endl;
69 if (vectorName == "Coordinates")
70 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Use CoordinatesTransferFactory to transfer coordinates instead of MultiVectorTransferFactory.");
71
72 RCP<MultiVector> fineVector = fineLevel.Get<RCP<MultiVector>>(vectorName, GetFactory("Vector factory").get());
73 RCP<MultiVector> coarseVector;
74
75 if (!isUncoupledAggFact) {
76 RCP<Matrix> transferOp = coarseLevel.Get<RCP<Matrix>>(transferName, GetFactory("Transfer factory").get());
77 Teuchos::ETransp transp;
78
79 if (transferOp->getGlobalNumRows() <= transferOp->getGlobalNumCols()) {
80 // R mode
81 coarseVector = MultiVectorFactory::Build(transferOp->getRangeMap(), fineVector->getNumVectors());
82 transp = Teuchos::NO_TRANS;
83 } else {
84 // P mode
85 coarseVector = MultiVectorFactory::Build(transferOp->getDomainMap(), fineVector->getNumVectors());
86 transp = Teuchos::TRANS;
87 }
88
89 transferOp->apply(*fineVector, *coarseVector, transp);
90
91 if (normalize) {
92 // Do constant row sum normalization
93 RCP<MultiVector> onesVector = MultiVectorFactory::Build(fineVector->getMap(), 1);
94 onesVector->putScalar(Teuchos::ScalarTraits<Scalar>::one());
95 RCP<MultiVector> rowSumVector = MultiVectorFactory::Build(coarseVector->getMap(), 1);
96 transferOp->apply(*onesVector, *rowSumVector, transp);
97
98 RCP<Vector> rowSumReciprocalVector = VectorFactory::Build(coarseVector->getMap(), 1);
99 rowSumReciprocalVector->reciprocal(*rowSumVector);
100
101 RCP<MultiVector> coarseVectorNormalized = MultiVectorFactory::Build(coarseVector->getMap(), fineVector->getNumVectors());
102 coarseVectorNormalized->elementWiseMultiply(1.0, *rowSumReciprocalVector, *coarseVector, 0.0);
103
104 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVectorNormalized);
105 } else {
106 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVector);
107 }
108 } else {
109 using execution_space = typename Node::execution_space;
110 using ATS = KokkosKernels::ArithTraits<Scalar>;
111 using impl_scalar_type = typename ATS::val_type;
112
113 auto aggregates = fineLevel.Get<RCP<Aggregates>>(transferName, GetFactory("Transfer factory").get());
114 TEUCHOS_ASSERT(!aggregates->AggregatesCrossProcessors());
115 RCP<const Map> coarseMap = Get<RCP<const Map>>(fineLevel, "CoarseMap");
116
117 auto aggGraph = aggregates->GetGraph();
118 auto numAggs = aggGraph.numRows();
119
120 RCP<const Map> coarseVectorMap;
121
122 LO blkSize = 1;
123 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
124 blkSize = rcp_dynamic_cast<const StridedMap>(coarseMap)->getFixedBlockSize();
125
126 if (blkSize == 1) {
127 // Scalar system
128 // No amalgamation required, we can use the coarseMap
129 coarseVectorMap = coarseMap;
130 } else {
131 // Vector system
132 AmalgamationFactory<SC, LO, GO, NO>::AmalgamateMap(rcp_dynamic_cast<const StridedMap>(coarseMap), coarseVectorMap);
133 }
134
135 coarseVector = MultiVectorFactory::Build(coarseVectorMap, fineVector->getNumVectors());
136
137 auto lcl_fineVector = fineVector->getLocalViewDevice(Tpetra::Access::ReadOnly);
138 auto lcl_coarseVector = coarseVector->getLocalViewDevice(Tpetra::Access::OverwriteAll);
139
140 Kokkos::parallel_for(
141 "MueLu:MultiVectorTransferFactory",
142 Kokkos::RangePolicy<LocalOrdinal, execution_space>(0, numAggs),
143 KOKKOS_LAMBDA(const LO i) {
144 auto aggregate = aggGraph.rowConst(i);
145 for (size_t j = 0; j < lcl_coarseVector.extent(1); j++) {
146 impl_scalar_type sum = 0.0;
147 for (size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
148 sum += lcl_fineVector(aggregate(colID), j);
149 if (normalize)
150 lcl_coarseVector(i, j) = sum / aggregate.length;
151 else
152 lcl_coarseVector(i, j) = sum;
153 }
154 });
155 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVector);
156 }
157} // Build
158
159} // namespace MueLu
160
161#endif // MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
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.
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.
Class that holds all level-specific information.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
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.
Namespace for MueLu classes and methods.
@ Runtime0
One-liner description of what is happening.