46#ifndef MUELU_REGIONRFACTORY_DEF_HPP
47#define MUELU_REGIONRFACTORY_DEF_HPP
49#include <Xpetra_Matrix.hpp>
50#include <Xpetra_CrsGraphFactory.hpp>
61 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 RCP<ParameterList> validParamList = rcp(
new ParameterList());
66 validParamList->set<RCP<const FactoryBase> >(
"A", Teuchos::null,
67 "Generating factory of the matrix A");
68 validParamList->set<RCP<const FactoryBase> >(
"numDimensions", Teuchos::null,
69 "Number of spatial dimensions in the problem.");
70 validParamList->set<RCP<const FactoryBase> >(
"lNodesPerDim", Teuchos::null,
71 "Number of local nodes per spatial dimension on the fine grid.");
72 validParamList->set<RCP<const FactoryBase> >(
"Nullspace", Teuchos::null,
73 "Fine level nullspace used to construct the coarse level nullspace.");
74 validParamList->set<RCP<const FactoryBase> >(
"Coordinates", Teuchos::null,
75 "Fine level coordinates used to construct piece-wise linear prolongator and coarse level coordinates.");
76 validParamList->set<
bool> (
"keep coarse coords",
false,
"Flag to keep coordinates for special coarse grid solve");
78 return validParamList;
81 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
85 Input(fineLevel,
"A");
86 Input(fineLevel,
"numDimensions");
87 Input(fineLevel,
"lNodesPerDim");
88 Input(fineLevel,
"Nullspace");
89 Input(fineLevel,
"Coordinates");
93 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
98 RCP<Teuchos::FancyOStream> out;
99 if(
const char* dbg = std::getenv(
"MUELU_REGIONRFACTORY_DEBUG")) {
100 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
101 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
103 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
106 *out <<
"Starting RegionRFactory::Build." << std::endl;
109 const int numDimensions = Get<int>(fineLevel,
"numDimensions");
110 Array<LO> lFineNodesPerDim(3, Teuchos::OrdinalTraits<LO>::one());
112 Array<LO> lNodesPerDim = Get<Array<LO> >(fineLevel,
"lNodesPerDim");
113 for(
int dim = 0; dim < numDimensions; ++dim) {
114 lFineNodesPerDim[dim] = lNodesPerDim[dim];
117 *out <<
"numDimensions " << numDimensions <<
" and lFineNodesPerDim: " << lFineNodesPerDim
121 if(numDimensions < 1 || numDimensions > 3) {
122 throw std::runtime_error(
"numDimensions must be 1, 2 or 3!");
124 for(
int dim = 0; dim < numDimensions; ++dim) {
125 if(lFineNodesPerDim[dim] % 3 != 1) {
126 throw std::runtime_error(
"The number of fine node in each direction need to be 3n+1");
129 Array<LO> lCoarseNodesPerDim(3, Teuchos::OrdinalTraits<LO>::one());
131 const RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel,
"A");
133 RCP<realvaluedmultivector_type> fineCoordinates, coarseCoordinates;
134 fineCoordinates = Get< RCP<realvaluedmultivector_type> >(fineLevel,
"Coordinates");
135 if(
static_cast<int>(fineCoordinates->getNumVectors()) != numDimensions) {
136 throw std::runtime_error(
"The number of vectors in the coordinates is not equal to numDimensions!");
144 if(numDimensions == 1) {
145 throw std::runtime_error(
"RegionRFactory no implemented for 1D case yet.");
146 }
else if(numDimensions == 2) {
147 throw std::runtime_error(
"RegionRFactory no implemented for 2D case yet.");
148 }
else if(numDimensions == 3) {
149 Build3D(numDimensions, lFineNodesPerDim, A, fineCoordinates,
150 R, coarseCoordinates, lCoarseNodesPerDim);
153 const Teuchos::ParameterList& pL = GetParameterList();
156 RCP<ParameterList> Tparams;
157 if(pL.isSublist(
"matrixmatrix: kernel params"))
158 Tparams=rcp(
new ParameterList(pL.sublist(
"matrixmatrix: kernel params")));
160 Tparams= rcp(
new ParameterList);
163 *out <<
"Compute P=R^t" << std::endl;
165 Tparams->set(
"compute global constants: temporaries",Tparams->get(
"compute global constants: temporaries",
false));
166 Tparams->set(
"compute global constants", Tparams->get(
"compute global constants",
false));
167 std::string label =
"MueLu::RegionR-transR" + Teuchos::toString(coarseLevel.
GetLevelID());
170 *out <<
"Compute coarse nullspace" << std::endl;
171 RCP<MultiVector> fineNullspace = Get<RCP<MultiVector> >(fineLevel,
"Nullspace");
172 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(R->getRowMap(),
173 fineNullspace->getNumVectors());
174 R->apply(*fineNullspace, *coarseNullspace, Teuchos::NO_TRANS, Teuchos::ScalarTraits<SC>::one(),
175 Teuchos::ScalarTraits<SC>::zero());
177 *out <<
"Set data on coarse level" << std::endl;
178 Set(coarseLevel,
"numDimensions", numDimensions);
179 Set(coarseLevel,
"lNodesPerDim", lCoarseNodesPerDim);
180 Set(coarseLevel,
"Nullspace", coarseNullspace);
181 Set(coarseLevel,
"Coordinates", coarseCoordinates);
182 if(pL.get<
bool>(
"keep coarse coords")) {
183 coarseLevel.
Set<RCP<realvaluedmultivector_type> >(
"Coordinates2", coarseCoordinates,
NoFactory::get());
186 R->SetFixedBlockSize(A->GetFixedBlockSize());
187 P->SetFixedBlockSize(A->GetFixedBlockSize());
189 Set(coarseLevel,
"R", R);
190 Set(coarseLevel,
"P", P);
194 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
196 Build3D(
const int numDimensions,
197 Teuchos::Array<LocalOrdinal>& lFineNodesPerDim,
198 const RCP<Matrix>& A,
202 Teuchos::Array<LocalOrdinal>& lCoarseNodesPerDim)
const {
203 using local_matrix_type =
typename CrsMatrix::local_matrix_type;
204 using local_graph_type =
typename local_matrix_type::staticcrsgraph_type;
205 using row_map_type =
typename local_matrix_type::row_map_type::non_const_type;
206 using entries_type =
typename local_matrix_type::index_type::non_const_type;
207 using values_type =
typename local_matrix_type::values_type::non_const_type;
210 RCP<Teuchos::FancyOStream> out;
211 if(
const char* dbg = std::getenv(
"MUELU_REGIONRFACTORY_DEBUG")) {
212 out = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
213 out->setShowAllFrontMatter(
false).setShowProcRank(
true);
215 out = Teuchos::getFancyOStream(rcp(
new Teuchos::oblackholestream()));
220 for(
int dim = 0; dim < numDimensions; ++dim) {
221 lCoarseNodesPerDim[dim] = lFineNodesPerDim[dim] / 3 + 1;
223 *out <<
"lCoarseNodesPerDim " << lCoarseNodesPerDim << std::endl;
226 const LO blkSize = A->GetFixedBlockSize();
227 *out <<
"blkSize " << blkSize << std::endl;
231 const LO numRows = blkSize*lCoarseNodesPerDim[0]*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[2];
232 const LO numCols = blkSize*lFineNodesPerDim[0]*lFineNodesPerDim[1]*lFineNodesPerDim[2];
237 RCP<Map> rowMap = MapFactory::Build(A->getRowMap()->lib(),
238 Teuchos::OrdinalTraits<GO>::invalid(),
240 A->getRowMap()->getIndexBase(),
241 A->getRowMap()->getComm());
243 RCP<Map> coordRowMap = MapFactory::Build(A->getRowMap()->lib(),
244 Teuchos::OrdinalTraits<GO>::invalid(),
245 lCoarseNodesPerDim[0]*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[2],
246 A->getRowMap()->getIndexBase(),
247 A->getRowMap()->getComm());
249 coarseCoordinates = Xpetra::MultiVectorFactory<real_type, LO, GO, NO>::Build(coordRowMap,
251 Array<ArrayRCP<const real_type> > fineCoordData(numDimensions);
252 Array<ArrayRCP<real_type> > coarseCoordData(numDimensions);
253 for(
int dim = 0; dim < numDimensions; ++dim) {
254 fineCoordData[dim] = fineCoordinates->getData(dim);
255 coarseCoordData[dim] = coarseCoordinates->getDataNonConst(dim);
262 const LO cornerStencilLength = 27;
263 const LO edgeStencilLength = 45;
264 const LO faceStencilLength = 75;
265 const LO interiorStencilLength = 125;
268 const LO numCorners = 8;
269 const LO numEdges = 4*(lCoarseNodesPerDim[0] - 2)
270 + 4*(lCoarseNodesPerDim[1] - 2)
271 + 4*(lCoarseNodesPerDim[2] - 2);
272 const LO numFaces = 2*(lCoarseNodesPerDim[0] - 2)*(lCoarseNodesPerDim[1] - 2)
273 + 2*(lCoarseNodesPerDim[0] - 2)*(lCoarseNodesPerDim[2] - 2)
274 + 2*(lCoarseNodesPerDim[1] - 2)*(lCoarseNodesPerDim[2] - 2);
275 const LO numInteriors = (lCoarseNodesPerDim[0] - 2)*(lCoarseNodesPerDim[1] - 2)
276 *(lCoarseNodesPerDim[2] - 2);
278 const LO nnz = (numCorners*cornerStencilLength + numEdges*edgeStencilLength
279 + numFaces*faceStencilLength + numInteriors*interiorStencilLength)*blkSize;
284 *out <<
"R statistics:" << std::endl
285 <<
" -numRows= " << numRows << std::endl
286 <<
" -numCols= " << numCols << std::endl
287 <<
" -nnz= " << nnz << std::endl;
289 row_map_type row_map(Kokkos::ViewAllocateWithoutInitializing(
"row_map"), numRows + 1);
290 typename row_map_type::HostMirror row_map_h = Kokkos::create_mirror_view(row_map);
292 entries_type entries(Kokkos::ViewAllocateWithoutInitializing(
"entries"), nnz);
293 typename entries_type::HostMirror entries_h = Kokkos::create_mirror_view(entries);
295 values_type values(Kokkos::ViewAllocateWithoutInitializing(
"values"), nnz);
296 typename values_type::HostMirror values_h = Kokkos::create_mirror_view(values);
301 Array<SC> coeffs({1.0/3.0, 2.0/3.0, 1.0, 2.0/3.0, 1.0/3.0});
306 const LO edgeLineOffset = 2*cornerStencilLength + (lCoarseNodesPerDim[0] - 2)*edgeStencilLength;
307 const LO faceLineOffset = 2*edgeStencilLength + (lCoarseNodesPerDim[0] - 2)*faceStencilLength;
308 const LO interiorLineOffset = 2*faceStencilLength
309 + (lCoarseNodesPerDim[0] - 2)*interiorStencilLength;
311 const LO facePlaneOffset = 2*edgeLineOffset + (lCoarseNodesPerDim[1] - 2)*faceLineOffset;
312 const LO interiorPlaneOffset = 2*faceLineOffset + (lCoarseNodesPerDim[1] - 2)*interiorLineOffset;
319 LO coordRowIdx = 0, rowIdx = 0, coordColumnOffset = 0, columnOffset = 0, entryOffset = 0;
320 for(LO l = 0; l < blkSize; ++l) {
321 for(LO k = 0; k < 3; ++k) {
322 for(LO j = 0; j < 3; ++j) {
323 for(LO i = 0; i < 3; ++i) {
324 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
325 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i)*blkSize + l;
326 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k + 2]*coeffs[j + 2]*coeffs[i + 2];
331 for(LO l = 0; l < blkSize; ++l) {
332 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
334 for(
int dim = 0; dim <numDimensions; ++dim) {
335 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
339 coordRowIdx += (lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0];
340 rowIdx = coordRowIdx*blkSize;
341 coordColumnOffset += (lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0];
342 columnOffset = coordColumnOffset*blkSize;
343 entryOffset += (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset)*blkSize;
344 for(LO l = 0; l < blkSize; ++l) {
345 for(LO k = 0; k < 3; ++k) {
346 for(LO j = 0; j < 3; ++j) {
347 for(LO i = 0; i < 3; ++i) {
348 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
349 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i)*blkSize + l;
350 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k]*coeffs[j + 2]*coeffs[i + 2];
355 for(LO l = 0; l < blkSize; ++l) {
356 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
358 for(
int dim = 0; dim <numDimensions; ++dim) {
359 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
363 coordRowIdx = (lCoarseNodesPerDim[0] - 1);
364 rowIdx = coordRowIdx*blkSize;
365 coordColumnOffset = (lFineNodesPerDim[0] - 1);
366 columnOffset = coordColumnOffset*blkSize;
367 entryOffset = (cornerStencilLength + (lCoarseNodesPerDim[0] - 2)*edgeStencilLength)*blkSize;
368 for(LO l = 0; l < blkSize; ++l) {
369 for(LO k = 0; k < 3; ++k) {
370 for(LO j = 0; j < 3; ++j) {
371 for(LO i = 0; i < 3; ++i) {
372 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
373 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0]
374 + j*lFineNodesPerDim[0] + (i - 2))*blkSize + l;
375 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k + 2]*coeffs[j + 2]*coeffs[i];
380 for(LO l = 0; l < blkSize; ++l) {
381 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
383 for(
int dim = 0; dim <numDimensions; ++dim) {
384 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
388 coordRowIdx += (lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0];
389 rowIdx = coordRowIdx*blkSize;
390 coordColumnOffset += (lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0];
391 columnOffset = coordColumnOffset*blkSize;
392 entryOffset += (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset)*blkSize;
393 for(LO l = 0; l < blkSize; ++l) {
394 for(LO k = 0; k < 3; ++k) {
395 for(LO j = 0; j < 3; ++j) {
396 for(LO i = 0; i < 3; ++i) {
397 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
398 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i - 2)*blkSize + l;
399 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k]*coeffs[j + 2]*coeffs[i];
404 for(LO l = 0; l < blkSize; ++l) {
405 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
407 for(
int dim = 0; dim <numDimensions; ++dim) {
408 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
412 coordRowIdx = (lCoarseNodesPerDim[1] - 1)*lCoarseNodesPerDim[0];
413 rowIdx = coordRowIdx*blkSize;
414 coordColumnOffset = (lFineNodesPerDim[1] - 1)*lFineNodesPerDim[0];
415 columnOffset = coordColumnOffset*blkSize;
416 entryOffset = (edgeLineOffset + (lCoarseNodesPerDim[1] - 2)*faceLineOffset)*blkSize;
417 for(LO l = 0; l < blkSize; ++l) {
418 for(LO k = 0; k < 3; ++k) {
419 for(LO j = 0; j < 3; ++j) {
420 for(LO i = 0; i < 3; ++i) {
421 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
422 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0]
423 + (j - 2)*lFineNodesPerDim[0] + i)*blkSize + l;
424 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k + 2]*coeffs[j]*coeffs[i + 2];
429 for(LO l = 0; l < blkSize; ++l) {
430 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
432 for(
int dim = 0; dim <numDimensions; ++dim) {
433 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
437 coordRowIdx += (lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0];
438 rowIdx = coordRowIdx*blkSize;
439 coordColumnOffset += (lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0];
440 columnOffset = coordColumnOffset*blkSize;
441 entryOffset += (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset)*blkSize;
442 for(LO l = 0; l < blkSize; ++l) {
443 for(LO k = 0; k < 3; ++k) {
444 for(LO j = 0; j < 3; ++j) {
445 for(LO i = 0; i < 3; ++i) {
446 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
447 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i)*blkSize + l;
448 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i + 2];
453 for(LO l = 0; l < blkSize; ++l) {
454 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
456 for(
int dim = 0; dim <numDimensions; ++dim) {
457 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
461 coordRowIdx = (lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0] - 1);
462 rowIdx = coordRowIdx*blkSize;
463 coordColumnOffset = (lFineNodesPerDim[1]*lFineNodesPerDim[0] - 1);
464 columnOffset = coordColumnOffset*blkSize;
465 entryOffset = (edgeLineOffset + (lCoarseNodesPerDim[1] - 2)*faceLineOffset +
466 cornerStencilLength + (lCoarseNodesPerDim[0] - 2)*edgeStencilLength)*blkSize;
467 for(LO l = 0; l < blkSize; ++l) {
468 for(LO k = 0; k < 3; ++k) {
469 for(LO j = 0; j < 3; ++j) {
470 for(LO i = 0; i < 3; ++i) {
471 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
472 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0]
473 + (j - 2)*lFineNodesPerDim[0] + (i - 2))*blkSize + l;
474 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k + 2]*coeffs[j]*coeffs[i];
479 for(LO l = 0; l < blkSize; ++l) {
480 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
482 for(
int dim = 0; dim <numDimensions; ++dim) {
483 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
487 coordRowIdx += (lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0];
488 rowIdx = coordRowIdx*blkSize;
489 coordColumnOffset += (lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0];
490 columnOffset = coordColumnOffset*blkSize;
491 entryOffset += (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset)*blkSize;
492 for(LO l = 0; l < blkSize; ++l) {
493 for(LO k = 0; k < 3; ++k) {
494 for(LO j = 0; j < 3; ++j) {
495 for(LO i = 0; i < 3; ++i) {
496 entries_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = columnOffset
497 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + (i - 2))*blkSize + l;
498 values_h(entryOffset + k*9 + j*3 + i + cornerStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
503 for(LO l = 0; l < blkSize; ++l) {
504 row_map_h(rowIdx + 1 + l) = entryOffset + cornerStencilLength*(l+1);
506 for(
int dim = 0; dim <numDimensions; ++dim) {
507 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
512 if(lCoarseNodesPerDim[0] - 2 > 0) {
514 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
515 for(LO edgeIdx = 0; edgeIdx < lCoarseNodesPerDim[0] - 2; ++edgeIdx) {
518 coordRowIdx = (edgeIdx + 1);
519 rowIdx = coordRowIdx*blkSize;
520 coordColumnOffset = (edgeIdx + 1)*3;
521 columnOffset = coordColumnOffset*blkSize;
522 entryOffset = (cornerStencilLength + edgeIdx*edgeStencilLength)*blkSize;
523 for(LO l = 0; l < blkSize; ++l) {
524 for(LO k = 0; k < 3; ++k) {
525 for(LO j = 0; j < 3; ++j) {
526 for(LO i = 0; i < 5; ++i) {
527 entries_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = columnOffset
528 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i - 2)*blkSize + l;
529 values_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = coeffs[k + 2]*coeffs[j + 2]*coeffs[i];
534 for(LO l = 0; l < blkSize; ++l) {
535 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
537 for(
int dim = 0; dim <numDimensions; ++dim) {
538 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
542 coordRowIdx = ((lCoarseNodesPerDim[1] - 1)*lCoarseNodesPerDim[0] + edgeIdx + 1);
543 rowIdx = coordRowIdx*blkSize;
544 coordColumnOffset = ((lFineNodesPerDim[1] - 1)*lFineNodesPerDim[0] + (edgeIdx + 1)*3);
545 columnOffset = coordColumnOffset*blkSize;
546 entryOffset = (edgeLineOffset + (lCoarseNodesPerDim[1] - 2)*faceLineOffset
547 + cornerStencilLength + edgeIdx*edgeStencilLength)*blkSize;
548 for(LO l = 0; l < blkSize; ++l) {
549 for(LO k = 0; k < 3; ++k) {
550 for(LO j = 0; j < 3; ++j) {
551 for(LO i = 0; i < 5; ++i) {
552 entries_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = columnOffset
553 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
554 values_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = coeffs[k + 2]*coeffs[j]*coeffs[i];
559 for(LO l = 0; l < blkSize; ++l) {
560 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
562 for(
int dim = 0; dim <numDimensions; ++dim) {
563 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
567 coordRowIdx = ((lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
569 rowIdx = coordRowIdx*blkSize;
570 coordColumnOffset = ((lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
572 columnOffset = coordColumnOffset*blkSize;
573 entryOffset = (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset
574 + cornerStencilLength + edgeIdx*edgeStencilLength)*blkSize;
575 for(LO l = 0; l < blkSize; ++l) {
576 for(LO k = 0; k < 3; ++k) {
577 for(LO j = 0; j < 3; ++j) {
578 for(LO i = 0; i < 5; ++i) {
579 entries_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = columnOffset
580 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i - 2)*blkSize + l;
581 values_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j + 2]*coeffs[i];
586 for(LO l = 0; l < blkSize; ++l) {
587 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
589 for(
int dim = 0; dim <numDimensions; ++dim) {
590 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
594 coordRowIdx = ((lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
595 + (lCoarseNodesPerDim[1] - 1)*lCoarseNodesPerDim[0] + edgeIdx + 1);
596 rowIdx = coordRowIdx*blkSize;
597 coordColumnOffset = ((lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
598 + (lFineNodesPerDim[1] - 1)*lFineNodesPerDim[0] + (edgeIdx + 1)*3);
599 columnOffset = coordColumnOffset*blkSize;
600 entryOffset = (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset
601 + edgeLineOffset + (lCoarseNodesPerDim[1] - 2)*faceLineOffset
602 + cornerStencilLength + edgeIdx*edgeStencilLength)*blkSize;
603 for(LO l = 0; l < blkSize; ++l) {
604 for(LO k = 0; k < 3; ++k) {
605 for(LO j = 0; j < 3; ++j) {
606 for(LO i = 0; i < 5; ++i) {
607 entries_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = columnOffset
608 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
609 + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
610 values_h(entryOffset + k*15 + j*5 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
615 for(LO l = 0; l < blkSize; ++l) {
616 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
618 for(
int dim = 0; dim <numDimensions; ++dim) {
619 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
625 if(lCoarseNodesPerDim[1] - 2 > 0) {
627 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
628 for(LO edgeIdx = 0; edgeIdx < lCoarseNodesPerDim[1] - 2; ++edgeIdx) {
631 coordRowIdx = (edgeIdx + 1)*lCoarseNodesPerDim[0];
632 rowIdx = coordRowIdx*blkSize;
633 coordColumnOffset = (edgeIdx + 1)*3*lFineNodesPerDim[0];
634 columnOffset = coordColumnOffset*blkSize;
635 entryOffset = (edgeLineOffset + edgeIdx*faceLineOffset)*blkSize;
636 for(LO l = 0; l < blkSize; ++l) {
637 for(LO k = 0; k < 3; ++k) {
638 for(LO j = 0; j < 5; ++j) {
639 for(LO i = 0; i < 3; ++i) {
640 entries_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = columnOffset
641 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i)*blkSize + l;
642 values_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = coeffs[k + 2]*coeffs[j]*coeffs[i + 2];
647 for(LO l = 0; l < blkSize; ++l) {
648 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
650 for(
int dim = 0; dim <numDimensions; ++dim) {
651 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
655 coordRowIdx = ((edgeIdx + 1)*lCoarseNodesPerDim[0] + lCoarseNodesPerDim[0] - 1);
656 rowIdx = coordRowIdx*blkSize;
657 coordColumnOffset = ((edgeIdx + 1)*3*lFineNodesPerDim[0] + lFineNodesPerDim[0] - 1);
658 columnOffset = coordColumnOffset*blkSize;
659 entryOffset = (edgeLineOffset + edgeIdx*faceLineOffset
660 + edgeStencilLength + (lCoarseNodesPerDim[0] - 2)*faceStencilLength)*blkSize;
661 for(LO l = 0; l < blkSize; ++l) {
662 for(LO k = 0; k < 3; ++k) {
663 for(LO j = 0; j < 5; ++j) {
664 for(LO i = 0; i < 3; ++i) {
665 entries_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = columnOffset
666 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
667 values_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = coeffs[k + 2]*coeffs[j]*coeffs[i];
672 for(LO l = 0; l < blkSize; ++l) {
673 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
675 for(
int dim = 0; dim <numDimensions; ++dim) {
676 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
680 coordRowIdx = ((lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
681 + (edgeIdx + 1)*lCoarseNodesPerDim[0]);
682 rowIdx = coordRowIdx*blkSize;
683 coordColumnOffset = ((lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
684 + (edgeIdx + 1)*3*lFineNodesPerDim[0]);
685 columnOffset = coordColumnOffset*blkSize;
686 entryOffset = (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset
687 + edgeLineOffset + edgeIdx*faceLineOffset)*blkSize;
688 for(LO l = 0; l < blkSize; ++l) {
689 for(LO k = 0; k < 3; ++k) {
690 for(LO j = 0; j < 5; ++j) {
691 for(LO i = 0; i < 3; ++i) {
692 entries_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = columnOffset
693 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i)*blkSize + l;
694 values_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i + 2];
699 for(LO l = 0; l < blkSize; ++l) {
700 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
702 for(
int dim = 0; dim <numDimensions; ++dim) {
703 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
707 coordRowIdx = ((lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
708 + (edgeIdx + 1)*lCoarseNodesPerDim[0] + lCoarseNodesPerDim[0] - 1);
709 rowIdx = coordRowIdx*blkSize;
710 coordColumnOffset = ((lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
711 + (edgeIdx + 1)*3*lFineNodesPerDim[0] + lFineNodesPerDim[0] - 1);
712 columnOffset = coordColumnOffset*blkSize;
713 entryOffset = (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset
714 + edgeLineOffset + edgeIdx*faceLineOffset
715 + edgeStencilLength + (lCoarseNodesPerDim[0] - 2)*faceStencilLength)*blkSize;
716 for(LO l = 0; l < blkSize; ++l) {
717 for(LO k = 0; k < 3; ++k) {
718 for(LO j = 0; j < 5; ++j) {
719 for(LO i = 0; i < 3; ++i) {
720 entries_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = columnOffset
721 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
722 + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
723 values_h(entryOffset + k*15 + j*3 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
728 for(LO l = 0; l < blkSize; ++l) {
729 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
731 for(
int dim = 0; dim <numDimensions; ++dim) {
732 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
738 if(lCoarseNodesPerDim[2] - 2 > 0) {
740 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
741 for(LO edgeIdx = 0; edgeIdx < lCoarseNodesPerDim[2] - 2; ++edgeIdx) {
744 coordRowIdx = (edgeIdx + 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0];
745 rowIdx = coordRowIdx*blkSize;
746 coordColumnOffset = (edgeIdx + 1)*3*lFineNodesPerDim[1]*lFineNodesPerDim[0];
747 columnOffset = coordColumnOffset*blkSize;
748 entryOffset = (facePlaneOffset + edgeIdx*interiorPlaneOffset)*blkSize;
749 for(LO l = 0; l < blkSize; ++l) {
750 for(LO k = 0; k < 5; ++k) {
751 for(LO j = 0; j < 3; ++j) {
752 for(LO i = 0; i < 3; ++i) {
753 entries_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = columnOffset
754 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i)*blkSize + l;
755 values_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j + 2]*coeffs[i + 2];
760 for(LO l = 0; l < blkSize; ++l) {
761 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
763 for(
int dim = 0; dim <numDimensions; ++dim) {
764 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
768 coordRowIdx = ((edgeIdx + 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
769 + lCoarseNodesPerDim[0] - 1);
770 rowIdx = coordRowIdx*blkSize;
771 coordColumnOffset = ((edgeIdx + 1)*3*lFineNodesPerDim[1]*lFineNodesPerDim[0]
772 + lFineNodesPerDim[0] - 1);
773 columnOffset = coordColumnOffset*blkSize;
774 entryOffset = (facePlaneOffset + faceLineOffset - edgeStencilLength
775 + edgeIdx*interiorPlaneOffset)*blkSize;
776 for(LO l = 0; l < blkSize; ++l) {
777 for(LO k = 0; k < 5; ++k) {
778 for(LO j = 0; j < 3; ++j) {
779 for(LO i = 0; i < 3; ++i) {
780 entries_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = columnOffset
781 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i - 2)*blkSize + l;
782 values_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j + 2]*coeffs[i];
787 for(LO l = 0; l < blkSize; ++l) {
788 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
790 for(
int dim = 0; dim <numDimensions; ++dim) {
791 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
795 coordRowIdx = ((edgeIdx + 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
796 + (lCoarseNodesPerDim[1] - 1)*lCoarseNodesPerDim[0]);
797 rowIdx = coordRowIdx*blkSize;
798 coordColumnOffset = ((edgeIdx + 1)*3*lFineNodesPerDim[1]*lFineNodesPerDim[0]
799 + (lFineNodesPerDim[1] - 1)*lFineNodesPerDim[0]);
800 columnOffset = coordColumnOffset*blkSize;
801 entryOffset = (facePlaneOffset + edgeIdx*interiorPlaneOffset + faceLineOffset
802 + (lCoarseNodesPerDim[1] - 2)*interiorLineOffset)*blkSize;
803 for(LO l = 0; l < blkSize; ++l) {
804 for(LO k = 0; k < 5; ++k) {
805 for(LO j = 0; j < 3; ++j) {
806 for(LO i = 0; i < 3; ++i) {
807 entries_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = columnOffset
808 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i)*blkSize + l;
809 values_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i + 2];
814 for(LO l = 0; l < blkSize; ++l) {
815 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
817 for(
int dim = 0; dim <numDimensions; ++dim) {
818 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
822 coordRowIdx = ((edgeIdx + 2)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0] - 1);
823 rowIdx = coordRowIdx*blkSize;
824 coordColumnOffset = ((edgeIdx + 1)*3*lFineNodesPerDim[1]*lFineNodesPerDim[0]
825 + lFineNodesPerDim[1]*lFineNodesPerDim[0] - 1);
826 columnOffset = coordColumnOffset*blkSize;
827 entryOffset = (facePlaneOffset + (edgeIdx + 1)*interiorPlaneOffset - edgeStencilLength)*blkSize;
828 for(LO l = 0; l < blkSize; ++l) {
829 for(LO k = 0; k < 5; ++k) {
830 for(LO j = 0; j < 3; ++j) {
831 for(LO i = 0; i < 3; ++i) {
832 entries_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = columnOffset
833 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
834 + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
835 values_h(entryOffset + k*9 + j*3 + i + edgeStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
840 for(LO l = 0; l < blkSize; ++l) {
841 row_map_h(rowIdx + 1 + l) = entryOffset + edgeStencilLength*(l+1);
843 for(
int dim = 0; dim <numDimensions; ++dim) {
844 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
850 if((lCoarseNodesPerDim[0] - 2 > 0) && (lCoarseNodesPerDim[1] - 2 > 0)) {
852 Array<LO> gridIdx(3);
853 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
854 for(LO faceIdx=0; faceIdx < (lCoarseNodesPerDim[1]-2)*(lCoarseNodesPerDim[0]-2); ++faceIdx) {
857 coordRowIdx = ((gridIdx[1] + 1)*lCoarseNodesPerDim[0] + gridIdx[0] + 1);
858 rowIdx = coordRowIdx*blkSize;
859 coordColumnOffset = 3*((gridIdx[1] + 1)*lFineNodesPerDim[0] + gridIdx[0] + 1);
860 columnOffset = coordColumnOffset*blkSize;
861 entryOffset = (edgeLineOffset + edgeStencilLength
862 + gridIdx[1]*faceLineOffset + gridIdx[0]*faceStencilLength)*blkSize;
863 for(LO l = 0; l < blkSize; ++l) {
864 for(LO k = 0; k < 3; ++k) {
865 for(LO j = 0; j < 5; ++j) {
866 for(LO i = 0; i < 5; ++i) {
867 entries_h(entryOffset + k*25 + j*5 + i + faceStencilLength*l) = columnOffset
868 + (k*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
869 values_h(entryOffset + k*25 + j*5 + i + faceStencilLength*l) = coeffs[k + 2]*coeffs[j]*coeffs[i];
874 for(LO l = 0; l < blkSize; ++l) {
875 row_map_h(rowIdx + 1 + l) = entryOffset + faceStencilLength*(l+1);
877 for(
int dim = 0; dim <numDimensions; ++dim) {
878 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
882 coordRowIdx += (lCoarseNodesPerDim[2] - 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0];
883 rowIdx = coordRowIdx*blkSize;
884 coordColumnOffset += (lFineNodesPerDim[2] - 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0];
885 columnOffset = coordColumnOffset*blkSize;
886 entryOffset += (facePlaneOffset + (lCoarseNodesPerDim[2] - 2)*interiorPlaneOffset)*blkSize;
887 for(LO l = 0; l < blkSize; ++l) {
888 for(LO k = 0; k < 3; ++k) {
889 for(LO j = 0; j < 5; ++j) {
890 for(LO i = 0; i < 5; ++i) {
891 entries_h(entryOffset + k*25 + j*5 + i + faceStencilLength*l) = columnOffset
892 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
893 + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
894 values_h(entryOffset + k*25 + j*5 + i + faceStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
899 for(LO l = 0; l < blkSize; ++l) {
900 row_map_h(rowIdx + 1 + l) = entryOffset + faceStencilLength*(l+1);
902 for(
int dim = 0; dim <numDimensions; ++dim) {
903 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
910 if(gridIdx[0] == lCoarseNodesPerDim[0] - 2) {
918 if((lCoarseNodesPerDim[0] - 2 > 0) && (lCoarseNodesPerDim[2] - 2 > 0)) {
920 Array<LO> gridIdx(3);
921 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
922 for(LO faceIdx=0; faceIdx < (lCoarseNodesPerDim[2]-2)*(lCoarseNodesPerDim[0]-2); ++faceIdx) {
925 coordRowIdx = ((gridIdx[2] + 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0] + (gridIdx[0] + 1));
926 rowIdx = coordRowIdx*blkSize;
927 coordColumnOffset = ((gridIdx[2] + 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
929 columnOffset = coordColumnOffset*blkSize;
930 entryOffset = (facePlaneOffset + gridIdx[2]*interiorPlaneOffset + edgeStencilLength
931 + gridIdx[0]*faceStencilLength)*blkSize;
932 for(LO l = 0; l < blkSize; ++l) {
933 for(LO k = 0; k < 5; ++k) {
934 for(LO j = 0; j < 3; ++j) {
935 for(LO i = 0; i < 5; ++i) {
936 entries_h(entryOffset + k*15 + j*5 + i + faceStencilLength*l) = columnOffset
937 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + j*lFineNodesPerDim[0] + i - 2)*blkSize + l;
938 values_h(entryOffset + k*15 + j*5 + i + faceStencilLength*l) = coeffs[k]*coeffs[j + 2]*coeffs[i];
943 for(LO l = 0; l < blkSize; ++l) {
944 row_map_h(rowIdx + 1 + l) = entryOffset + faceStencilLength*(l+1);
946 for(
int dim = 0; dim <numDimensions; ++dim) {
947 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
951 coordRowIdx += (lCoarseNodesPerDim[1] - 1)*lCoarseNodesPerDim[0];
952 rowIdx = coordRowIdx*blkSize;
953 coordColumnOffset += (lFineNodesPerDim[1] - 1)*lFineNodesPerDim[0];
954 columnOffset = coordColumnOffset*blkSize;
955 entryOffset += (faceLineOffset + (lCoarseNodesPerDim[1] - 2)*interiorLineOffset)*blkSize;
956 for(LO l = 0; l < blkSize; ++l) {
957 for(LO k = 0; k < 5; ++k) {
958 for(LO j = 0; j < 3; ++j) {
959 for(LO i = 0; i < 5; ++i) {
960 entries_h(entryOffset + k*15 + j*5 + i + faceStencilLength*l) = columnOffset
961 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
962 + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
963 values_h(entryOffset + k*15 + j*5 + i + faceStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
968 for(LO l = 0; l < blkSize; ++l) {
969 row_map_h(rowIdx + 1 + l) = entryOffset + faceStencilLength*(l+1);
971 for(
int dim = 0; dim <numDimensions; ++dim) {
972 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
979 if(gridIdx[0] == lCoarseNodesPerDim[0] - 2) {
987 if((lCoarseNodesPerDim[1] - 2 > 0) && (lCoarseNodesPerDim[2] - 2 > 0)) {
989 Array<LO> gridIdx(3);
990 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
991 for(LO faceIdx=0; faceIdx < (lCoarseNodesPerDim[2]-2)*(lCoarseNodesPerDim[1]-2); ++faceIdx) {
994 coordRowIdx = ((gridIdx[2] + 1)*lCoarseNodesPerDim[1]*lCoarseNodesPerDim[0]
995 + (gridIdx[1] + 1)*lCoarseNodesPerDim[0]);
996 rowIdx = coordRowIdx*blkSize;
997 coordColumnOffset = ((gridIdx[2] + 1)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
998 + (gridIdx[1] + 1)*lFineNodesPerDim[0])*3;
999 columnOffset = coordColumnOffset*blkSize;
1000 entryOffset = (facePlaneOffset + gridIdx[2]*interiorPlaneOffset + faceLineOffset
1001 + gridIdx[1]*interiorLineOffset)*blkSize;
1002 for(LO l = 0; l < blkSize; ++l) {
1003 for(LO k = 0; k < 5; ++k) {
1004 for(LO j = 0; j < 5; ++j) {
1005 for(LO i = 0; i < 3; ++i) {
1006 entries_h(entryOffset + k*15 + j*3 + i + faceStencilLength*l) = columnOffset
1007 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0] + (j - 2)*lFineNodesPerDim[0] + i)*blkSize + l;
1008 values_h(entryOffset + k*15 + j*3 + i + faceStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i + 2];
1013 for(LO l = 0; l < blkSize; ++l) {
1014 row_map_h(rowIdx + 1 + l) = entryOffset + faceStencilLength*(l+1);
1016 for(
int dim = 0; dim <numDimensions; ++dim) {
1017 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
1021 coordRowIdx += (lCoarseNodesPerDim[0] - 1);
1022 rowIdx = coordRowIdx*blkSize;
1023 coordColumnOffset += (lFineNodesPerDim[0] - 1);
1024 columnOffset = coordColumnOffset*blkSize;
1025 entryOffset += (faceStencilLength + (lCoarseNodesPerDim[0] - 2)*interiorStencilLength)*blkSize;
1026 for(LO l = 0; l < blkSize; ++l) {
1027 for(LO k = 0; k < 5; ++k) {
1028 for(LO j = 0; j < 5; ++j) {
1029 for(LO i = 0; i < 3; ++i) {
1030 entries_h(entryOffset + k*15 + j*3 + i + faceStencilLength*l) = columnOffset
1031 + ((k - 2)*lFineNodesPerDim[1]*lFineNodesPerDim[0]
1032 + (j - 2)*lFineNodesPerDim[0] + i - 2)*blkSize + l;
1033 values_h(entryOffset + k*15 + j*3 + i + faceStencilLength*l) = coeffs[k]*coeffs[j]*coeffs[i];
1038 for(LO l = 0; l < blkSize; ++l) {
1039 row_map_h(rowIdx + 1 + l) = entryOffset + faceStencilLength*(l+1);
1041 for(
int dim = 0; dim <numDimensions; ++dim) {
1042 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
1049 if(gridIdx[1] == lCoarseNodesPerDim[1] - 2) {
1056 if(numInteriors > 0) {
1061 LO countRowEntries = 0;
1062 Array<LO> coordColumnOffsets(125);
1063 for(LO k = -2; k < 3; ++k) {
1064 for(LO j = -2; j < 3; ++j) {
1065 for(LO i = -2; i < 3; ++i) {
1066 coordColumnOffsets[countRowEntries] = k*lFineNodesPerDim[1]*lFineNodesPerDim[0]
1067 + j*lFineNodesPerDim[0] + i;
1074 Array<SC> interiorValues(125);
1075 for(LO k = 0; k < 5; ++k) {
1076 for(LO j = 0; j < 5; ++j) {
1077 for(LO i = 0; i < 5; ++i) {
1078 interiorValues[countValues] = coeffs[k]*coeffs[j]*coeffs[i];
1084 LO coordRowIdx, rowIdx, coordColumnOffset, columnOffset, entryOffset;
1085 Array<LO> gridIdx(3);
1086 for(LO interiorIdx = 0; interiorIdx < numInteriors; ++interiorIdx) {
1087 coordRowIdx = ((gridIdx[2] + 1)*lCoarseNodesPerDim[0]*lCoarseNodesPerDim[1]
1088 + (gridIdx[1] + 1)*lCoarseNodesPerDim[0]
1090 rowIdx = coordRowIdx*blkSize;
1091 coordColumnOffset = ((gridIdx[2] + 1)*3*lFineNodesPerDim[1]*lFineNodesPerDim[0]
1092 + (gridIdx[1] + 1)*3*lFineNodesPerDim[0] + (gridIdx[0] + 1)*3);
1093 columnOffset = coordColumnOffset*blkSize;
1094 entryOffset = (facePlaneOffset + faceLineOffset + faceStencilLength
1095 + gridIdx[2]*interiorPlaneOffset + gridIdx[1]*interiorLineOffset
1096 + gridIdx[0]*interiorStencilLength)*blkSize;
1097 for(LO l = 0; l < blkSize; ++l) {
1098 row_map_h(rowIdx + 1 + l) = entryOffset + interiorStencilLength*(l+1);
1103 for(LO l = 0; l < blkSize; ++l) {
1104 for(LO entryIdx = 0; entryIdx < interiorStencilLength; ++entryIdx) {
1105 entries_h(entryOffset + entryIdx + interiorStencilLength*l) = columnOffset + coordColumnOffsets[entryIdx]*blkSize + l;
1106 values_h(entryOffset + entryIdx + interiorStencilLength*l) = interiorValues[entryIdx];
1109 for(
int dim = 0; dim <numDimensions; ++dim) {
1110 coarseCoordData[dim][coordRowIdx] = fineCoordData[dim][coordColumnOffset];
1117 if(gridIdx[0] == lCoarseNodesPerDim[0] - 2) {
1120 if(gridIdx[1] == lCoarseNodesPerDim[1] - 2) {
1128 Kokkos::deep_copy(row_map, row_map_h);
1129 Kokkos::deep_copy(entries, entries_h);
1130 Kokkos::deep_copy(values, values_h);
1132 local_graph_type localGraph(entries, row_map);
1133 local_matrix_type localR(
"R", numCols, values, localGraph);
1135 R = MatrixFactory::Build(localR,
1146#define MUELU_REGIONRFACTORY_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
static const NoFactory * get()
void Build3D(const int numDimensions, Array< LO > &lFineNodesPerDim, const RCP< Matrix > &A, const RCP< realvaluedmultivector_type > &fineCoordinates, RCP< Matrix > &R, RCP< realvaluedmultivector_type > &coarseCoordinates, Array< LO > &lCoarseNodesPerDim) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
RCP< const ParameterList > GetValidParameterList() const
Input.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Namespace for MueLu classes and methods.