MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_CombinePFactory_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_COMBINEPFACTORY_DEF_HPP
11#define MUELU_COMBINEPFACTORY_DEF_HPP
12
13#include <stdlib.h>
14#include <iomanip>
15
16// #include <Teuchos_LAPACK.hpp>
17#include <Teuchos_SerialDenseMatrix.hpp>
18#include <Teuchos_SerialDenseVector.hpp>
19#include <Teuchos_SerialDenseSolver.hpp>
20
21#include <Xpetra_CrsMatrixWrap.hpp>
22#include <Xpetra_ImportFactory.hpp>
23#include <Xpetra_Matrix.hpp>
24#include <Xpetra_MapFactory.hpp>
25#include <Xpetra_MultiVectorFactory.hpp>
26#include <Xpetra_VectorFactory.hpp>
27
28#include <Xpetra_IO.hpp>
29
32
33#include "MueLu_MasterList.hpp"
34#include "MueLu_Monitor.hpp"
35
36namespace MueLu {
37
38template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
40 RCP<ParameterList> validParamList = rcp(new ParameterList());
41 validParamList->setEntry("combine: numBlks", ParameterEntry(1));
42 validParamList->set<RCP<const FactoryBase>>("A", Teuchos::null, "Generating factory of the matrix A");
43
44 return validParamList;
45}
46
47template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
49 DeclareInput(Level& fineLevel, Level& /* coarseLevel */) const {
50 // Input(fineLevel, "subblock");
51}
52
53template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
55 Level& coarseLevel) const {
56 return BuildP(fineLevel, coarseLevel);
57}
58
59template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61 Level& coarseLevel) const {
62 FactoryMonitor m(*this, "Build", coarseLevel);
63
64 const ParameterList& pL = GetParameterList();
65 const LO nBlks = as<LO>(pL.get<int>("combine: numBlks"));
66
67 RCP<Matrix> A = Get<RCP<Matrix>>(fineLevel, "A");
68
69 // Record all matrices that each define a block in block diagonal comboP
70 // matrix used for PDE/multiblock interpolation. Additionally, count
71 // total number of local rows, nonzeros, coarseDofs, and colDofs.
72
73 Teuchos::ArrayRCP<RCP<Matrix>> arrayOfMatrices(nBlks);
74 size_t nComboRowMap = 0, nnzCombo = 0, nComboColMap = 0, nComboDomMap = 0;
75
76 LO nTotalNumberLocalColMapEntries = 0;
77 Teuchos::ArrayRCP<size_t> DomMapSizePerBlk(nBlks);
78 Teuchos::ArrayRCP<size_t> ColMapSizePerBlk(nBlks);
79 Teuchos::ArrayRCP<size_t> ColMapLocalSizePerBlk(nBlks);
80 Teuchos::ArrayRCP<size_t> ColMapRemoteSizePerBlk(nBlks);
81 Teuchos::ArrayRCP<size_t> ColMapLocalCumulativePerBlk(nBlks + 1); // hardwire 0th entry so that it has the value of 0
82 Teuchos::ArrayRCP<size_t> ColMapRemoteCumulativePerBlk(nBlks + 1); // hardwire 0th entry so that it has the value of 0
83 for (int j = 0; j < nBlks; j++) {
84 std::string blockName = "Psubblock" + Teuchos::toString(j);
85 if (coarseLevel.IsAvailable(blockName, NoFactory::get())) {
86 arrayOfMatrices[j] = coarseLevel.Get<RCP<Matrix>>(blockName, NoFactory::get());
87 nComboRowMap += Teuchos::as<size_t>((arrayOfMatrices[j])->getRowMap()->getLocalNumElements());
88 DomMapSizePerBlk[j] = Teuchos::as<size_t>((arrayOfMatrices[j])->getDomainMap()->getLocalNumElements());
89 ColMapSizePerBlk[j] = Teuchos::as<size_t>((arrayOfMatrices[j])->getColMap()->getLocalNumElements());
90 nComboDomMap += DomMapSizePerBlk[j];
91 nComboColMap += ColMapSizePerBlk[j];
92 nnzCombo += Teuchos::as<size_t>((arrayOfMatrices[j])->getLocalNumEntries());
93 TEUCHOS_TEST_FOR_EXCEPTION((arrayOfMatrices[j])->getDomainMap()->getIndexBase() != 0, Exceptions::RuntimeError, "interpolation subblocks must use 0 indexbase");
94
95 // figure out how many empty entries in each column map
96 int tempii = 0;
97 for (int i = 0; i < (int)DomMapSizePerBlk[j]; i++) {
98 // if ( (arrayOfMatrices[j])->getDomainMap()->getGlobalElement(i) == (arrayOfMatrices[j])->getColMap()->getGlobalElement(tempii) ) nTotalNumberLocalColMapEntries++;
99 if ((arrayOfMatrices[j])->getDomainMap()->getGlobalElement(i) == (arrayOfMatrices[j])->getColMap()->getGlobalElement(tempii)) tempii++;
100 }
101 nTotalNumberLocalColMapEntries += tempii;
102 ColMapLocalSizePerBlk[j] = tempii;
103 ColMapRemoteSizePerBlk[j] = ColMapSizePerBlk[j] - ColMapLocalSizePerBlk[j];
104 } else {
105 arrayOfMatrices[j] = Teuchos::null;
106 ColMapLocalSizePerBlk[j] = 0;
107 ColMapRemoteSizePerBlk[j] = 0;
108 }
109 ColMapLocalCumulativePerBlk[j + 1] = ColMapLocalSizePerBlk[j] + ColMapLocalCumulativePerBlk[j];
110 ColMapRemoteCumulativePerBlk[j + 1] = ColMapRemoteSizePerBlk[j] + ColMapRemoteCumulativePerBlk[j];
111 }
112 TEUCHOS_TEST_FOR_EXCEPTION(nComboRowMap != A->getRowMap()->getLocalNumElements(), Exceptions::RuntimeError, "sum of subblock rows != #row's Afine");
113
114 // build up csr arrays for combo block diagonal P
115 Teuchos::ArrayRCP<size_t> comboPRowPtr(nComboRowMap + 1);
116 Teuchos::ArrayRCP<LocalOrdinal> comboPCols(nnzCombo);
117 Teuchos::ArrayRCP<Scalar> comboPVals(nnzCombo);
118
119 size_t nnzCnt = 0, nrowCntFromPrevBlks = 0, ncolCntFromPrevBlks = 0;
120 LO maxNzPerRow = 0;
121 for (int j = 0; j < nBlks; j++) {
122 // grab csr pointers for individual blocks of P
123 if (arrayOfMatrices[j] != Teuchos::null) {
124 Teuchos::ArrayRCP<const size_t> subblockRowPtr((arrayOfMatrices[j])->getLocalNumRows());
125 Teuchos::ArrayRCP<const LocalOrdinal> subblockCols((arrayOfMatrices[j])->getLocalNumEntries());
126 Teuchos::ArrayRCP<const Scalar> subblockVals((arrayOfMatrices[j])->getLocalNumEntries());
127 if ((int)(arrayOfMatrices[j])->getLocalMaxNumRowEntries() > maxNzPerRow) maxNzPerRow = (int)(arrayOfMatrices[j])->getLocalMaxNumRowEntries();
128 Teuchos::RCP<CrsMatrixWrap> subblockwrap = Teuchos::rcp_dynamic_cast<CrsMatrixWrap>((arrayOfMatrices[j]));
129 Teuchos::RCP<CrsMatrix> subblockcrs = subblockwrap->getCrsMatrix();
130 subblockcrs->getAllValues(subblockRowPtr, subblockCols, subblockVals);
131
132 // copy jth block into csr arrays of comboP
133
134 for (decltype(subblockRowPtr.size()) i = 0; i < subblockRowPtr.size() - 1; i++) {
135 size_t rowLength = subblockRowPtr[i + 1] - subblockRowPtr[i];
136 comboPRowPtr[nrowCntFromPrevBlks + i] = nnzCnt;
137 for (size_t k = 0; k < rowLength; k++) {
138 if ((int)subblockCols[k + subblockRowPtr[i]] < (int)ColMapLocalSizePerBlk[j]) {
139 comboPCols[nnzCnt] = subblockCols[k + subblockRowPtr[i]] + ColMapLocalCumulativePerBlk[j];
140 if ((int)comboPCols[nnzCnt] >= (int)ColMapLocalCumulativePerBlk[nBlks]) {
141 printf("ERROR1\n");
142 exit(1);
143 }
144 } else {
145 // Here we subtract off the number of local colmap guys ... so this tell us where we are among ghost unknowns. We then want to stick this ghost after
146 // all the Local guys ... so we add ColMapLocalCumulativePerBlk[nBlks] .... finally we need to account for the fact that other ghost blocks may have already
147 // been handled ... so we then add + ColMapRemoteCumulativePerBlk[j];
148 comboPCols[nnzCnt] = subblockCols[k + subblockRowPtr[i]] - ColMapLocalSizePerBlk[j] + ColMapLocalCumulativePerBlk[nBlks] + ColMapRemoteCumulativePerBlk[j];
149 if ((int)comboPCols[nnzCnt] >= (int)(ColMapLocalCumulativePerBlk[nBlks] + ColMapRemoteCumulativePerBlk[nBlks])) {
150 printf("ERROR2\n");
151 exit(1);
152 }
153 }
154 comboPVals[nnzCnt++] = subblockVals[k + subblockRowPtr[i]];
155 }
156 }
157
158 nrowCntFromPrevBlks += Teuchos::as<size_t>((arrayOfMatrices[j])->getRowMap()->getLocalNumElements());
159 ncolCntFromPrevBlks += DomMapSizePerBlk[j]; // rst: check this
160 }
161 }
162 comboPRowPtr[nComboRowMap] = nnzCnt;
163
164 // Come up with global IDs for the coarse grid maps. We assume that each xxx
165 // block has a minimum GID of 0. Since MueLu is generally creating these
166 // GIDS, this assumption is probably correct, but we'll check it.
167
168 Teuchos::Array<GlobalOrdinal> comboDomainMapGIDs(nComboDomMap);
169 Teuchos::Array<GlobalOrdinal> comboColMapGIDs(nComboColMap);
170
171 LO nTotalNumberRemoteColMapEntries = 0;
172 GlobalOrdinal offset = 0;
173 size_t domainMapIndex = 0;
174 int nComboColIndex = 0;
175 for (int j = 0; j < nBlks; j++) {
176 int nThisBlkColIndex = 0;
177
178 GlobalOrdinal tempMax = 0, maxGID = 0;
179 if (arrayOfMatrices[j] != Teuchos::null) tempMax = (arrayOfMatrices[j])->getDomainMap()->getMaxGlobalIndex();
180 Teuchos::reduceAll(*(A->getDomainMap()->getComm()), Teuchos::REDUCE_MAX, tempMax, Teuchos::ptr(&maxGID));
181
182 if (arrayOfMatrices[j] != Teuchos::null) {
183 TEUCHOS_TEST_FOR_EXCEPTION(arrayOfMatrices[j]->getDomainMap()->getMinAllGlobalIndex() < 0, Exceptions::RuntimeError, "Global ID assumption for domainMap not met within subblock");
184
185 GO priorDomGID = 0;
186 for (size_t c = 0; c < DomMapSizePerBlk[j]; ++c) { // check this
187 // The global ids of jth block are assumed to go between 0 and maxGID_j. So the 1th blocks DomainGIDs should start at maxGID_0+1. The 2nd
188 // block DomainDIGS starts at maxGID_0+1 + maxGID_1 + 1. We use offset to keep track of these starting points.
189 priorDomGID = (arrayOfMatrices[j])->getDomainMap()->getGlobalElement(c);
190 comboDomainMapGIDs[domainMapIndex++] = offset + priorDomGID;
191 if (priorDomGID == (arrayOfMatrices[j])->getColMap()->getGlobalElement(nThisBlkColIndex)) {
192 comboColMapGIDs[nComboColIndex++] = offset + priorDomGID;
193 nThisBlkColIndex++;
194 }
195 }
196
197 for (size_t cc = nThisBlkColIndex; cc < ColMapSizePerBlk[j]; ++cc) {
198 comboColMapGIDs[nTotalNumberLocalColMapEntries + nTotalNumberRemoteColMapEntries++] = offset + (arrayOfMatrices[j])->getColMap()->getGlobalElement(cc);
199 }
200 }
201 offset += maxGID + 1;
202 }
203#ifdef out
204 RCP<const Map> coarseDomainMap = Xpetra::MapFactory<LO, GO, NO>::Build(A->getDomainMap()->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), comboDomainMapGIDs, 0, A->getDomainMap()->getComm());
205 RCP<const Map> coarseColMap = Xpetra::MapFactory<LO, GO, NO>::Build(A->getDomainMap()->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), comboColMapGIDs, 0, A->getDomainMap()->getComm());
206#endif
207
208 RCP<const Map> coarseDomainMap = Xpetra::MapFactory<LO, GO, NO>::Build(A->getDomainMap()->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), comboDomainMapGIDs, 0, A->getRowMap()->getComm());
209 RCP<const Map> coarseColMap = Xpetra::MapFactory<LO, GO, NO>::Build(A->getDomainMap()->lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), comboColMapGIDs, 0, A->getRowMap()->getComm());
210
211 Teuchos::RCP<CrsMatrix> comboPCrs = CrsMatrixFactory::Build(A->getRowMap(), coarseColMap, maxNzPerRow);
212 // comboPCrs->getCrsGraph(); //.getRowInfo(6122);
213 // comboPCrs->getRowInfo(6122);
214
215 // Teuchos::RCP<CrsMatrix> comboPCrs = CrsMatrixFactory::Build(A->getRowMap(), coarseColMap,nnzCombo+1000);
216
217 // for (size_t i = 0; i < nComboRowMap; i++) {
218 // printf("FIXME\n"); if (nComboRowMap > 6142) nComboRowMap = 6142;
219 for (size_t i = 0; i < nComboRowMap; i++) {
220 comboPCrs->insertLocalValues(i, comboPCols.view(comboPRowPtr[i], comboPRowPtr[i + 1] - comboPRowPtr[i]),
221 comboPVals.view(comboPRowPtr[i], comboPRowPtr[i + 1] - comboPRowPtr[i]));
222 }
223 comboPCrs->fillComplete(coarseDomainMap, A->getRowMap());
224
225 Teuchos::RCP<Matrix> comboP = Teuchos::rcp(new CrsMatrixWrap(comboPCrs));
226
227 Set(coarseLevel, "P", comboP);
228}
229
230} // namespace MueLu
231
232#define MUELU_COMBINEPFACTORY_SHORT
233#endif // MUELU_COMBINEPFACTORY_DEF_HPP
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
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.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
static const NoFactory * get()
Namespace for MueLu classes and methods.