10 #ifndef MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
11 #define MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
14 #include "Xpetra_Access.hpp"
15 #include "Xpetra_MultiVectorFactory.hpp"
18 #include "MueLu_UncoupledAggregationFactory.hpp"
19 #include "MueLu_Aggregates.hpp"
24 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
28 validParamList->
set<std::string>(
"Vector name",
"undefined",
"Name of the vector that will be transferred on the coarse grid (level key)");
29 validParamList->
set<std::string>(
"Transfer name",
"R",
"Name of the operator that will be used to transfer data");
30 validParamList->
set<
bool>(
"Normalize",
false,
"If a row sum normalization should be applied to preserve the mean value of the vector.");
35 return validParamList;
38 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
41 std::string vectorName = pL.
get<std::string>(
"Vector name");
42 std::string transferName = pL.
get<std::string>(
"Transfer name");
44 fineLevel.
DeclareInput(vectorName, GetFactory(
"Vector factory").
get(),
this);
45 auto transferFact = GetFactory(
"Transfer factory");
47 if (isUncoupledAggFact) {
48 fineLevel.
DeclareInput(transferName, transferFact.get(),
this);
49 Input(fineLevel,
"CoarseMap");
51 coarseLevel.
DeclareInput(transferName, transferFact.get(),
this);
54 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
59 const std::string vectorName = pL.
get<std::string>(
"Vector name");
60 const std::string transferName = pL.
get<std::string>(
"Transfer name");
61 const bool normalize = pL.
get<
bool>(
"Normalize");
63 auto transferFact = GetFactory(
"Transfer factory");
66 GetOStream(
Runtime0) <<
"Transferring multivector \"" << vectorName <<
"\" using operator \"" << transferName <<
"\"" << std::endl;
67 if (vectorName ==
"Coordinates")
73 if (!isUncoupledAggFact) {
77 if (transferOp->getGlobalNumRows() <= transferOp->getGlobalNumCols()) {
79 coarseVector = MultiVectorFactory::Build(transferOp->getRangeMap(), fineVector->getNumVectors());
83 coarseVector = MultiVectorFactory::Build(transferOp->getDomainMap(), fineVector->getNumVectors());
87 transferOp->apply(*fineVector, *coarseVector, transp);
91 RCP<MultiVector> onesVector = MultiVectorFactory::Build(fineVector->getMap(), 1);
93 RCP<MultiVector> rowSumVector = MultiVectorFactory::Build(coarseVector->getMap(), 1);
94 transferOp->apply(*onesVector, *rowSumVector, transp);
96 RCP<Vector> rowSumReciprocalVector = VectorFactory::Build(coarseVector->getMap(), 1);
97 rowSumReciprocalVector->reciprocal(*rowSumVector);
99 RCP<MultiVector> coarseVectorNormalized = MultiVectorFactory::Build(coarseVector->getMap(), fineVector->getNumVectors());
100 coarseVectorNormalized->elementWiseMultiply(1.0, *rowSumReciprocalVector, *coarseVector, 0.0);
102 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVectorNormalized);
104 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVector);
107 using execution_space =
typename Node::execution_space;
108 using ATS = Kokkos::ArithTraits<Scalar>;
109 using impl_scalar_type =
typename ATS::val_type;
111 auto aggregates = fineLevel.
Get<
RCP<Aggregates>>(transferName, GetFactory(
"Transfer factory").get());
113 RCP<const Map> coarseMap = Get<RCP<const Map>>(fineLevel,
"CoarseMap");
115 auto aggGraph = aggregates->GetGraph();
116 auto numAggs = aggGraph.numRows();
118 coarseVector = MultiVectorFactory::Build(coarseMap, fineVector->getNumVectors());
120 auto lcl_fineVector = fineVector->getDeviceLocalView(Xpetra::Access::ReadOnly);
121 auto lcl_coarseVector = coarseVector->getDeviceLocalView(Xpetra::Access::OverwriteAll);
123 Kokkos::parallel_for(
124 "MueLu:MultiVectorTransferFactory",
125 Kokkos::RangePolicy<LocalOrdinal, execution_space>(0, numAggs),
126 KOKKOS_LAMBDA(
const LO i) {
127 auto aggregate = aggGraph.rowConst(i);
128 for (
size_t j = 0; j < lcl_coarseVector.extent(1); j++) {
129 impl_scalar_type sum = 0.0;
130 for (
size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
131 sum += lcl_fineVector(aggregate(colID), j);
133 lcl_coarseVector(i, j) = sum / aggregate.length;
135 lcl_coarseVector(i, j) = sum;
138 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVector);
144 #endif // MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level->Get< RCP<Matrix> >("A", factory) if factory == NULL => use default factory.
bool is_null(const boost::shared_ptr< T > &p)
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Exception throws to report errors in the internal logical of the program.
#define TEUCHOS_ASSERT(assertion_test)
Factory for building uncoupled aggregates.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.