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 "Xpetra_Access.hpp"
15#include "Xpetra_MultiVectorFactory.hpp"
16
17#include "MueLu_Level.hpp"
18#include "MueLu_UncoupledAggregationFactory.hpp"
19#include "MueLu_Aggregates.hpp"
20#include "MueLu_Monitor.hpp"
21
22namespace MueLu {
23
24template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
26 RCP<ParameterList> validParamList = rcp(new ParameterList());
27
28 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?
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.");
31 validParamList->set<RCP<const FactoryBase>>("Vector factory", Teuchos::null, "Factory of the vector");
32 validParamList->set<RCP<const FactoryBase>>("Transfer factory", Teuchos::null, "Factory of the transfer operator");
33 validParamList->set<RCP<const FactoryBase>>("CoarseMap", Teuchos::null, "Generating factory of the coarse map");
34
35 return validParamList;
36}
37
38template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40 const ParameterList &pL = GetParameterList();
41 std::string vectorName = pL.get<std::string>("Vector name");
42 std::string transferName = pL.get<std::string>("Transfer name");
43
44 fineLevel.DeclareInput(vectorName, GetFactory("Vector factory").get(), this);
45 auto transferFact = GetFactory("Transfer factory");
46 const bool isUncoupledAggFact = !Teuchos::rcp_dynamic_cast<const UncoupledAggregationFactory>(transferFact).is_null();
47 if (isUncoupledAggFact) {
48 fineLevel.DeclareInput(transferName, transferFact.get(), this);
49 Input(fineLevel, "CoarseMap");
50 } else
51 coarseLevel.DeclareInput(transferName, transferFact.get(), this);
52}
53
54template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56 FactoryMonitor m(*this, "Build", coarseLevel);
57
58 const ParameterList &pL = GetParameterList();
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");
62
63 auto transferFact = GetFactory("Transfer factory");
64 const bool isUncoupledAggFact = !Teuchos::rcp_dynamic_cast<const UncoupledAggregationFactory>(transferFact).is_null();
65
66 GetOStream(Runtime0) << "Transferring multivector \"" << vectorName << "\" using operator \"" << transferName << "\"" << std::endl;
67 if (vectorName == "Coordinates")
68 TEUCHOS_TEST_FOR_EXCEPTION(true, Exceptions::RuntimeError, "Use CoordinatesTransferFactory to transfer coordinates instead of MultiVectorTransferFactory.");
69
70 RCP<MultiVector> fineVector = fineLevel.Get<RCP<MultiVector>>(vectorName, GetFactory("Vector factory").get());
71 RCP<MultiVector> coarseVector;
72
73 if (!isUncoupledAggFact) {
74 RCP<Matrix> transferOp = coarseLevel.Get<RCP<Matrix>>(transferName, GetFactory("Transfer factory").get());
75 Teuchos::ETransp transp;
76
77 if (transferOp->getGlobalNumRows() <= transferOp->getGlobalNumCols()) {
78 // R mode
79 coarseVector = MultiVectorFactory::Build(transferOp->getRangeMap(), fineVector->getNumVectors());
80 transp = Teuchos::NO_TRANS;
81 } else {
82 // P mode
83 coarseVector = MultiVectorFactory::Build(transferOp->getDomainMap(), fineVector->getNumVectors());
84 transp = Teuchos::TRANS;
85 }
86
87 transferOp->apply(*fineVector, *coarseVector, transp);
88
89 if (normalize) {
90 // Do constant row sum normalization
91 RCP<MultiVector> onesVector = MultiVectorFactory::Build(fineVector->getMap(), 1);
92 onesVector->putScalar(Teuchos::ScalarTraits<Scalar>::one());
93 RCP<MultiVector> rowSumVector = MultiVectorFactory::Build(coarseVector->getMap(), 1);
94 transferOp->apply(*onesVector, *rowSumVector, transp);
95
96 RCP<Vector> rowSumReciprocalVector = VectorFactory::Build(coarseVector->getMap(), 1);
97 rowSumReciprocalVector->reciprocal(*rowSumVector);
98
99 RCP<MultiVector> coarseVectorNormalized = MultiVectorFactory::Build(coarseVector->getMap(), fineVector->getNumVectors());
100 coarseVectorNormalized->elementWiseMultiply(1.0, *rowSumReciprocalVector, *coarseVector, 0.0);
101
102 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVectorNormalized);
103 } else {
104 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVector);
105 }
106 } else {
107 using execution_space = typename Node::execution_space;
108 using ATS = Kokkos::ArithTraits<Scalar>;
109 using impl_scalar_type = typename ATS::val_type;
110
111 auto aggregates = fineLevel.Get<RCP<Aggregates>>(transferName, GetFactory("Transfer factory").get());
112 TEUCHOS_ASSERT(!aggregates->AggregatesCrossProcessors());
113 RCP<const Map> coarseMap = Get<RCP<const Map>>(fineLevel, "CoarseMap");
114
115 auto aggGraph = aggregates->GetGraph();
116 auto numAggs = aggGraph.numRows();
117
118 coarseVector = MultiVectorFactory::Build(coarseMap, fineVector->getNumVectors());
119
120 auto lcl_fineVector = fineVector->getDeviceLocalView(Xpetra::Access::ReadOnly);
121 auto lcl_coarseVector = coarseVector->getDeviceLocalView(Xpetra::Access::OverwriteAll);
122
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);
132 if (normalize)
133 lcl_coarseVector(i, j) = sum / aggregate.length;
134 else
135 lcl_coarseVector(i, j) = sum;
136 }
137 });
138 Set<RCP<MultiVector>>(coarseLevel, vectorName, coarseVector);
139 }
140} // Build
141
142} // namespace MueLu
143
144#endif // MUELU_MULTIVECTORTRANSFER_FACTORY_DEF_HPP
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.