MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_UnsmooshFactory_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 PACKAGES_MUELU_SRC_GRAPH_MUELU_UNSMOOSHFACTORY_DEF_HPP_
11 #define PACKAGES_MUELU_SRC_GRAPH_MUELU_UNSMOOSHFACTORY_DEF_HPP_
12 
13 #include "MueLu_Monitor.hpp"
14 
16 
17 namespace MueLu {
18 
19 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
21 
22 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
24  RCP<ParameterList> validParamList = rcp(new ParameterList());
25  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory for unamalgamated matrix. Row map of (unamalgamted) output prolongation operator should match row map of this A.");
26  validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Generating factory of the (amalgamated) prolongator P");
27  validParamList->set<RCP<const FactoryBase> >("DofStatus", Teuchos::null, "Generating factory for dofStatus array (usually the VariableDofLaplacdianFactory)");
28 
29  validParamList->set<int>("maxDofPerNode", 1, "Maximum number of DOFs per node");
30  validParamList->set<bool>("fineIsPadded", false, "true if finest level input matrix is padded");
31 
32  return validParamList;
33 }
34 
35 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37  // const ParameterList& pL = GetParameterList();
38  Input(fineLevel, "A");
39  Input(coarseLevel, "P");
40 
41  // DofStatus only provided on the finest level (by user)
42  // On the coarser levels it is auto-generated using the DBC information from the unamalgamated matrix A
43  if (fineLevel.GetLevelID() == 0)
44  Input(fineLevel, "DofStatus");
45 }
46 
47 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
49  FactoryMonitor m(*this, "Build", coarseLevel);
50  typedef Teuchos::ScalarTraits<SC> STS;
51 
52  const ParameterList &pL = GetParameterList();
53 
54  // extract matrices (unamalgamated A and amalgamated P)
55  RCP<Matrix> unamalgA = Get<RCP<Matrix> >(fineLevel, "A");
56  RCP<Matrix> amalgP = Get<RCP<Matrix> >(coarseLevel, "P");
57 
58  // extract user parameters
59  int maxDofPerNode = pL.get<int>("maxDofPerNode");
60  bool fineIsPadded = pL.get<bool>("fineIsPadded");
61 
62  // get dofStatus information
63  // On the finest level it is provided by the user. On the coarser levels it is constructed
64  // using the DBC information of the matrix A
65  Teuchos::Array<char> dofStatus;
66  if (fineLevel.GetLevelID() == 0) {
67  dofStatus = Get<Teuchos::Array<char> >(fineLevel, "DofStatus");
68  } else {
69  // dof status is the dirichlet information of unsmooshed/unamalgamated A (fine level)
70  dofStatus = Teuchos::Array<char>(unamalgA->getRowMap()->getLocalNumElements() /*amalgP->getRowMap()->getLocalNumElements() * maxDofPerNode*/, 's');
71 
72  bool bHasZeroDiagonal = false;
74 
75  TEUCHOS_TEST_FOR_EXCEPTION(dirOrNot.size() != dofStatus.size(), MueLu::Exceptions::RuntimeError, "MueLu::UnsmooshFactory::Build: inconsistent number of coarse DBC array and dofStatus array. dirOrNot.size() = " << dirOrNot.size() << " dofStatus.size() = " << dofStatus.size());
76  for (decltype(dirOrNot.size()) i = 0; i < dirOrNot.size(); ++i) {
77  if (dirOrNot[i] == true) dofStatus[i] = 'p';
78  }
79  }
80 
81  // TODO: TAW the following check is invalid for SA-AMG based input prolongators
82  // TEUCHOS_TEST_FOR_EXCEPTION(amalgP->getDomainMap()->isSameAs(*amalgP->getColMap()) == false, MueLu::Exceptions::RuntimeError,"MueLu::UnsmooshFactory::Build: only support for non-overlapping aggregates. (column map of Ptent must be the same as domain map of Ptent)");
83 
84  // extract CRS information from amalgamated prolongation operator
85  Teuchos::ArrayRCP<const size_t> amalgRowPtr(amalgP->getLocalNumRows());
86  Teuchos::ArrayRCP<const LocalOrdinal> amalgCols(amalgP->getLocalNumEntries());
87  Teuchos::ArrayRCP<const Scalar> amalgVals(amalgP->getLocalNumEntries());
88  Teuchos::RCP<CrsMatrix> amalgPcrs = toCrsMatrix(amalgP);
89  amalgPcrs->getAllValues(amalgRowPtr, amalgCols, amalgVals);
90 
91  // calculate number of dof rows for new prolongator
92  size_t paddedNrows = amalgP->getRowMap()->getLocalNumElements() * Teuchos::as<size_t>(maxDofPerNode);
93 
94  // reserve CSR arrays for new prolongation operator
95  Teuchos::ArrayRCP<size_t> newPRowPtr(paddedNrows + 1);
96  Teuchos::ArrayRCP<LocalOrdinal> newPCols(amalgP->getLocalNumEntries() * maxDofPerNode);
97  Teuchos::ArrayRCP<Scalar> newPVals(amalgP->getLocalNumEntries() * maxDofPerNode);
98 
99  size_t rowCount = 0; // actual number of (local) in unamalgamated prolongator
100  if (fineIsPadded == true || fineLevel.GetLevelID() > 0) {
101  // build prolongation operator for padded fine level matrices.
102  // Note: padded fine level dofs are transferred by injection.
103  // That is, these interpolation stencils do not take averages of
104  // coarse level variables. Further, fine level Dirichlet points
105  // also use injection.
106 
107  size_t cnt = 0; // local id counter
108  for (decltype(amalgRowPtr.size()) i = 0; i < amalgRowPtr.size() - 1; i++) {
109  // determine number of entries in amalgamated dof row i
110  size_t rowLength = amalgRowPtr[i + 1] - amalgRowPtr[i];
111 
112  // loop over dofs per node (unamalgamation)
113  for (int j = 0; j < maxDofPerNode; j++) {
114  newPRowPtr[i * maxDofPerNode + j] = cnt;
115  if (dofStatus[i * maxDofPerNode + j] == 's') { // add only "standard" dofs to unamalgamated prolongator
116  // loop over column entries in amalgamated P
117  for (size_t k = 0; k < rowLength; k++) {
118  newPCols[cnt] = amalgCols[k + amalgRowPtr[i]] * maxDofPerNode + j;
119  newPVals[cnt++] = amalgVals[k + amalgRowPtr[i]];
120  }
121  }
122  }
123  }
124 
125  newPRowPtr[paddedNrows] = cnt; // close row CSR array
126  rowCount = paddedNrows;
127  } else {
128  // Build prolongation operator for non-padded fine level matrices.
129  // Need to map from non-padded dofs to padded dofs. For this, look
130  // at the status array and skip padded dofs.
131 
132  size_t cnt = 0; // local id counter
133 
134  for (decltype(amalgRowPtr.size()) i = 0; i < amalgRowPtr.size() - 1; i++) {
135  // determine number of entries in amalgamated dof row i
136  size_t rowLength = amalgRowPtr[i + 1] - amalgRowPtr[i];
137 
138  // loop over dofs per node (unamalgamation)
139  for (int j = 0; j < maxDofPerNode; j++) {
140  // no interpolation for padded fine dofs as they do not exist
141 
142  if (dofStatus[i * maxDofPerNode + j] == 's') { // add only "standard" dofs to unamalgamated prolongator
143  newPRowPtr[rowCount++] = cnt;
144  // loop over column entries in amalgamated P
145  for (size_t k = 0; k < rowLength; k++) {
146  newPCols[cnt] = amalgCols[k + amalgRowPtr[i]] * maxDofPerNode + j;
147  newPVals[cnt++] = amalgVals[k + amalgRowPtr[i]];
148  }
149  }
150  if (dofStatus[i * maxDofPerNode + j] == 'd') { // Dirichlet handling
151  newPRowPtr[rowCount++] = cnt;
152  }
153  }
154  }
155  newPRowPtr[rowCount] = cnt; // close row CSR array
156  } // fineIsPadded == false
157 
158  // generate coarse domain map
159  // So far no support for gid offset or strided maps. This information
160  // could be gathered easily from the unamalgamated fine level operator A.
161  std::vector<size_t> stridingInfo(1, maxDofPerNode);
162 
163  GlobalOrdinal nCoarseDofs = amalgP->getDomainMap()->getLocalNumElements() * maxDofPerNode;
164  GlobalOrdinal indexBase = amalgP->getDomainMap()->getIndexBase();
165  RCP<const Map> coarseDomainMap = StridedMapFactory::Build(amalgP->getDomainMap()->lib(),
167  nCoarseDofs,
168  indexBase,
169  stridingInfo,
170  amalgP->getDomainMap()->getComm(),
171  -1 /* stridedBlockId */,
172  0 /*domainGidOffset */);
173 
174  size_t nColCoarseDofs = Teuchos::as<size_t>(amalgP->getColMap()->getLocalNumElements() * maxDofPerNode);
175  Teuchos::Array<GlobalOrdinal> unsmooshColMapGIDs(nColCoarseDofs);
176  for (size_t c = 0; c < amalgP->getColMap()->getLocalNumElements(); ++c) {
177  GlobalOrdinal gid = (amalgP->getColMap()->getGlobalElement(c) - indexBase) * maxDofPerNode + indexBase;
178 
179  for (int i = 0; i < maxDofPerNode; ++i) {
180  unsmooshColMapGIDs[c * maxDofPerNode + i] = gid + i;
181  }
182  }
183  Teuchos::RCP<Map> coarseColMap = MapFactory::Build(amalgP->getDomainMap()->lib(),
185  unsmooshColMapGIDs(), // View,
186  indexBase,
187  amalgP->getDomainMap()->getComm());
188 
189  // Assemble unamalgamated P
190  Teuchos::RCP<CrsMatrix> unamalgPCrs = CrsMatrixFactory::Build(unamalgA->getRowMap(),
191  coarseColMap,
192  maxDofPerNode * amalgP->getLocalMaxNumRowEntries());
193  for (size_t i = 0; i < rowCount; i++) {
194  unamalgPCrs->insertLocalValues(i,
195  newPCols.view(newPRowPtr[i], newPRowPtr[i + 1] - newPRowPtr[i]),
196  newPVals.view(newPRowPtr[i], newPRowPtr[i + 1] - newPRowPtr[i]));
197  }
198  unamalgPCrs->fillComplete(coarseDomainMap, unamalgA->getRowMap());
199 
200  Teuchos::RCP<Matrix> unamalgP = Teuchos::rcp(new CrsMatrixWrap(unamalgPCrs));
201 
202  Set(coarseLevel, "P", unamalgP);
203 }
204 
205 } // namespace MueLu
206 
207 #endif /* PACKAGES_MUELU_SRC_GRAPH_MUELU_UNSMOOSHFACTORY_DEF_HPP_ */
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)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
size_type size() const
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
static Teuchos::ArrayRCP< const bool > DetectDirichletRowsExt(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool &bHasZeroDiagonal, const Magnitude &tol=Teuchos::ScalarTraits< Scalar >::zero())
Detect Dirichlet rows (extended version)
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:51
Exception throws to report errors in the internal logical of the program.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.