46 #ifndef MUELU_GEOINTERPFACTORY_DEF_HPP
47 #define MUELU_GEOINTERPFACTORY_DEF_HPP
54 #include <Xpetra_Map.hpp>
55 #include <Xpetra_MapFactory.hpp>
56 #include <Xpetra_Matrix.hpp>
57 #include <Xpetra_MultiVector.hpp>
58 #include <Xpetra_MultiVectorFactory.hpp>
59 #include <Xpetra_VectorFactory.hpp>
60 #include <Xpetra_CrsMatrixWrap.hpp>
61 #include <Xpetra_MultiVector.hpp>
62 #include <Xpetra_MultiVectorFactory.hpp>
67 #include <MueLu_Utilities.hpp>
73 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
75 GetOStream(
Runtime1) <<
"I constructed a GeoInterpFactory object... Nothing else to do here." << std::endl;
78 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
83 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
86 Input(fineLevel,
"A");
87 Input(fineLevel,
"A00");
88 Input(fineLevel,
"A10");
89 Input(fineLevel,
"A20");
91 Input(fineLevel,
"VElementList");
92 Input(fineLevel,
"PElementList");
93 Input(fineLevel,
"MElementList");
96 Input(coarseLevel,
"VElementList");
97 Input(coarseLevel,
"PElementList");
98 Input(coarseLevel,
"MElementList");
113 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
116 GetOStream(
Runtime1) <<
"Starting 'build' routine...\n";
124 return BuildP(fineLevel,coarseLevel);
128 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
133 GetOStream(
Runtime1) <<
"Starting 'BuildP' routine...\n";
145 GetOStream(
Runtime1) <<
"done getting fine level elements...\n";
146 GetOStream(
Runtime1) <<
"getting coarse level elements...\n";
155 coarseLevel.
Get(
"VElementList",coarseElementVDOFs,coarseLevel.
GetFactoryManager()->GetFactory(
"VElementList").get());
156 coarseLevel.
Get(
"PElementList",coarseElementPDOFs,coarseLevel.
GetFactoryManager()->GetFactory(
"PElementList").get());
157 coarseLevel.
Get(
"MElementList",coarseElementMDOFs,coarseLevel.
GetFactoryManager()->GetFactory(
"MElementList").get());
160 GetOStream(
Runtime1) <<
"computing various numbers...\n";
162 GO totalFineElements = fineElementMDOFs->numRows();
163 LO nFineElements = (int) sqrt(totalFineElements);
164 GO totalCoarseElements = coarseElementMDOFs->numRows();
165 LO nCoarseElements = (int) sqrt(totalCoarseElements);
168 GO nM = (2*nCoarseElements+1)*(2*nCoarseElements+1);
170 GO nP = (nCoarseElements+1)*(nCoarseElements+1);
173 RCP<Matrix> fineA00 = Get<RCP<Matrix> > (fineLevel,
"A00");
174 RCP<Matrix> fineA10 = Get<RCP<Matrix> > (fineLevel,
"A10");
175 RCP<Matrix> fineA20 = Get<RCP<Matrix> > (fineLevel,
"A20");
177 GetOStream(
Runtime1) <<
"creating coarse grid maps...\n";
183 GO fNV = rowMapforPV->getGlobalNumElements();
184 GO fNP = rowMapforPP->getGlobalNumElements();
185 GO fNM = rowMapforPM->getGlobalNumElements();
191 RCP<Matrix> FineA = Factory::Get< RCP<Matrix> >(fineLevel,
"A");
195 RCP<const Map> colMapforPV = Xpetra::MapFactory<LO,GO>::createUniformContigMap(Xpetra::UseTpetra,nV,comm);
196 RCP<const Map> colMapforPP = Xpetra::MapFactory<LO,GO>::createUniformContigMap(Xpetra::UseTpetra,nP,comm);
197 RCP<const Map> colMapforPM = Xpetra::MapFactory<LO,GO>::createUniformContigMap(Xpetra::UseTpetra,nM,comm);
200 GetOStream(
Runtime1) <<
"creating coarse grid matrices...\n";
202 size_t maxEntriesPerRowV = 9,
203 maxEntriesPerRowP = 4,
204 maxEntriesPerRowM = 9;
206 RCP<Matrix> P =
rcp(
new CrsMatrixWrap(rowMapforP ,maxEntriesPerRowV));
207 RCP<Matrix> PV =
rcp(
new CrsMatrixWrap(rowMapforPV,maxEntriesPerRowV));
208 RCP<Matrix> PP =
rcp(
new CrsMatrixWrap(rowMapforPP,maxEntriesPerRowP));
209 RCP<Matrix> PM =
rcp(
new CrsMatrixWrap(rowMapforPM,maxEntriesPerRowM));
269 GO fineElement[4] = {0,1,nFineElements,nFineElements+1};
271 GetOStream(
Runtime1) <<
"start building matrices...\n";
272 GetOStream(
Runtime1) <<
"nCoarseElements = " << nCoarseElements << std::endl;
274 for ( GO coarseElement=0; coarseElement<totalCoarseElements; coarseElement++)
284 if (coarseElement < nCoarseElements)
288 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
290 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,1);
292 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,4);
296 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],4),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
298 colPtrV[0] = 2*colPtrM[0];
299 colPtrV[1] = 2*colPtrM[1];
300 colPtrV[2] = 2*colPtrM[2];
302 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],8),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
304 colPtrV[0] = 2*colPtrM[0]+1;
305 colPtrV[1] = 2*colPtrM[1]+1;
306 colPtrV[2] = 2*colPtrM[2]+1;
308 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],9),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
311 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,0);
313 colPtrP[1] = (*coarseElementPDOFs)(coarseElement,1);
317 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[0],1),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
320 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,4);
324 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],1),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
326 colPtrV[0] = 2*colPtrM[0];
328 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],2),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
330 colPtrV[0] = 2*colPtrM[0]+1;
332 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],3),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
335 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,1);
339 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[1],1),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
343 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
345 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,1);
347 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,4);
351 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[1],4),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
353 colPtrV[0] = 2*colPtrM[0];
354 colPtrV[1] = 2*colPtrM[1];
355 colPtrV[2] = 2*colPtrM[2];
357 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],8),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
359 colPtrV[0] = 2*colPtrM[0]+1;
360 colPtrV[1] = 2*colPtrM[1]+1;
361 colPtrV[2] = 2*colPtrM[2]+1;
363 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],9),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
367 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,1);
371 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[1],1),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
373 colPtrV[0] = 2*colPtrM[0];
375 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],2),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
377 colPtrV[0] = 2*colPtrM[0]+1;
379 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],3),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
382 if (coarseElement == 0)
387 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
391 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],0),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
393 colPtrV[0] = 2*colPtrM[0];
395 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],0),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
397 colPtrV[0] = 2*colPtrM[0]+1;
399 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],1),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
402 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,0);
406 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[0],0),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
412 if (coarseElement % (nCoarseElements) == 0)
416 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
418 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,3);
420 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,7);
424 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],7),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
426 colPtrV[0] = 2*colPtrM[0];
427 colPtrV[1] = 2*colPtrM[1];
428 colPtrV[2] = 2*colPtrM[2];
430 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],14),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
432 colPtrV[0] = 2*colPtrM[0]+1;
433 colPtrV[1] = 2*colPtrM[1]+1;
434 colPtrV[2] = 2*colPtrM[2]+1;
436 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],15),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
439 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,7);
443 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],3),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
445 colPtrV[0] = 2*colPtrM[0];
447 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],6),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
449 colPtrV[0] = 2*colPtrM[0]+1;
451 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],7),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
454 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
456 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,3);
458 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,7);
462 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[2],7),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
464 colPtrV[0] = 2*colPtrM[0];
465 colPtrV[1] = 2*colPtrM[1];
466 colPtrV[2] = 2*colPtrM[2];
468 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],14),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
470 colPtrV[0] = 2*colPtrM[0]+1;
471 colPtrV[1] = 2*colPtrM[1]+1;
472 colPtrV[2] = 2*colPtrM[2]+1;
474 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],15),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
478 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,3);
482 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[2],3),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
484 colPtrV[0] = 2*colPtrM[0];
486 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],6),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
488 colPtrV[0] = 2*colPtrM[0]+1;
490 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],7),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
495 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,0);
497 colPtrP[1] = (*coarseElementPDOFs)(coarseElement,3);
501 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[0],3),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
504 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,3);
508 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[2],3),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
515 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,8);
519 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],0),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
521 colPtrV[0] = 2*colPtrM[0];
523 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],0),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
525 colPtrV[0] = 2*colPtrM[0]+1;
527 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],1),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
530 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,5);
534 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],1),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
536 colPtrV[0] = 2*colPtrM[0];
538 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],2),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
540 colPtrV[0] = 2*colPtrM[0]+1;
542 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],3),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
545 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,2);
549 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],2),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
551 colPtrV[0] = 2*colPtrM[0];
553 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],4),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
555 colPtrV[0] = 2*colPtrM[0]+1;
557 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],5),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
560 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,6);
564 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],3),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
566 colPtrV[0] = 2*colPtrM[0];
568 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],6),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
570 colPtrV[0] = 2*colPtrM[0]+1;
572 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],7),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
575 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,4);
577 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,6);
579 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,8);
583 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],5),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
585 colPtrV[0] = 2*colPtrM[0];
586 colPtrV[1] = 2*colPtrM[1];
587 colPtrV[2] = 2*colPtrM[2];
589 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],10),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
591 colPtrV[0] = 2*colPtrM[0]+1;
592 colPtrV[1] = 2*colPtrM[1]+1;
593 colPtrV[2] = 2*colPtrM[2]+1;
595 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],11),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
599 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,5);
601 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,7);
603 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,8);
607 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],6),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
609 colPtrV[0] = 2*colPtrM[0];
610 colPtrV[1] = 2*colPtrM[1];
611 colPtrV[2] = 2*colPtrM[2];
613 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],12),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
615 colPtrV[0] = 2*colPtrM[0]+1;
616 colPtrV[1] = 2*colPtrM[1]+1;
617 colPtrV[2] = 2*colPtrM[2]+1;
619 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],13),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
623 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,1);
625 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,2);
627 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,5);
631 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[1],5),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
633 colPtrV[0] = 2*colPtrM[0];
634 colPtrV[1] = 2*colPtrM[1];
635 colPtrV[2] = 2*colPtrM[2];
637 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],10),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
639 colPtrV[0] = 2*colPtrM[0]+1;
640 colPtrV[1] = 2*colPtrM[1]+1;
641 colPtrV[2] = 2*colPtrM[2]+1;
643 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],11),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
647 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,5);
649 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,7);
651 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,8);
655 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[1],6),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
657 colPtrV[0] = 2*colPtrM[0];
658 colPtrV[1] = 2*colPtrM[1];
659 colPtrV[2] = 2*colPtrM[2];
661 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],12),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
663 colPtrV[0] = 2*colPtrM[0]+1;
664 colPtrV[1] = 2*colPtrM[1]+1;
665 colPtrV[2] = 2*colPtrM[2]+1;
667 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],13),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
671 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,4);
673 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,6);
675 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,8);
679 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[2],5),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
681 colPtrV[0] = 2*colPtrM[0];
682 colPtrV[1] = 2*colPtrM[1];
683 colPtrV[2] = 2*colPtrM[2];
685 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],10),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
687 colPtrV[0] = 2*colPtrM[0]+1;
688 colPtrV[1] = 2*colPtrM[1]+1;
689 colPtrV[2] = 2*colPtrM[2]+1;
691 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],11),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
695 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,2);
697 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,3);
699 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,6);
703 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[2],6),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
705 colPtrV[0] = 2*colPtrM[0];
706 colPtrV[1] = 2*colPtrM[1];
707 colPtrV[2] = 2*colPtrM[2];
709 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],12),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
711 colPtrV[0] = 2*colPtrM[0]+1;
712 colPtrV[1] = 2*colPtrM[1]+1;
713 colPtrV[2] = 2*colPtrM[2]+1;
715 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],13),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
719 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,1);
721 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,2);
723 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,5);
727 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],5),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
729 colPtrV[0] = 2*colPtrM[0];
730 colPtrV[1] = 2*colPtrM[1];
731 colPtrV[2] = 2*colPtrM[2];
733 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],10),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
735 colPtrV[0] = 2*colPtrM[0]+1;
736 colPtrV[1] = 2*colPtrM[1]+1;
737 colPtrV[2] = 2*colPtrM[2]+1;
740 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],11),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
743 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,2);
745 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,3);
747 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,6);
751 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],6),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
753 colPtrV[0] = 2*colPtrM[0];
754 colPtrV[1] = 2*colPtrM[1];
755 colPtrV[2] = 2*colPtrM[2];
757 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],12),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
759 colPtrV[0] = 2*colPtrM[0]+1;
760 colPtrV[1] = 2*colPtrM[1]+1;
761 colPtrV[2] = 2*colPtrM[2]+1;
763 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],13),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
767 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
768 valPtrM[0] = 0.140625;
769 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,1);
770 valPtrM[1] = -0.046875;
771 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,2);
772 valPtrM[2] = 0.015625;
773 colPtrM[3] = (*coarseElementMDOFs)(coarseElement,3);
774 valPtrM[3] = -0.046875;
775 colPtrM[4] = (*coarseElementMDOFs)(coarseElement,4);
776 valPtrM[4] = 0.28125;
777 colPtrM[5] = (*coarseElementMDOFs)(coarseElement,5);
778 valPtrM[5] = -0.09375;
779 colPtrM[6] = (*coarseElementMDOFs)(coarseElement,6);
780 valPtrM[6] = -0.09375;
781 colPtrM[7] = (*coarseElementMDOFs)(coarseElement,7);
782 valPtrM[7] = 0.28125;
783 colPtrM[8] = (*coarseElementMDOFs)(coarseElement,8);
787 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[0],8),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
789 colPtrV[0] = 2*colPtrM[0];
790 colPtrV[1] = 2*colPtrM[1];
791 colPtrV[2] = 2*colPtrM[2];
792 colPtrV[3] = 2*colPtrM[3];
793 colPtrV[4] = 2*colPtrM[4];
794 colPtrV[5] = 2*colPtrM[5];
795 colPtrV[6] = 2*colPtrM[6];
796 colPtrV[7] = 2*colPtrM[7];
797 colPtrV[8] = 2*colPtrM[8];
799 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],16),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
801 colPtrV[0] = 2*colPtrM[0]+1;
802 colPtrV[1] = 2*colPtrM[1]+1;
803 colPtrV[2] = 2*colPtrM[2]+1;
804 colPtrV[3] = 2*colPtrM[3]+1;
805 colPtrV[4] = 2*colPtrM[4]+1;
806 colPtrV[5] = 2*colPtrM[5]+1;
807 colPtrV[6] = 2*colPtrM[6]+1;
808 colPtrV[7] = 2*colPtrM[7]+1;
809 colPtrV[8] = 2*colPtrM[8]+1;
811 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[0],17),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
814 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
815 valPtrM[0] = -0.046875;
816 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,1);
817 valPtrM[1] = 0.140625;
818 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,2);
819 valPtrM[2] = -0.046875;
820 colPtrM[3] = (*coarseElementMDOFs)(coarseElement,3);
821 valPtrM[3] = 0.015625;
822 colPtrM[4] = (*coarseElementMDOFs)(coarseElement,4);
823 valPtrM[4] = 0.28125;
824 colPtrM[5] = (*coarseElementMDOFs)(coarseElement,5);
825 valPtrM[5] = 0.28125;
826 colPtrM[6] = (*coarseElementMDOFs)(coarseElement,6);
827 valPtrM[6] = -0.09375;
828 colPtrM[7] = (*coarseElementMDOFs)(coarseElement,7);
829 valPtrM[7] = -0.09375;
830 colPtrM[8] = (*coarseElementMDOFs)(coarseElement,8);
834 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[1],8),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
836 colPtrV[0] = 2*colPtrM[0];
837 colPtrV[1] = 2*colPtrM[1];
838 colPtrV[2] = 2*colPtrM[2];
839 colPtrV[3] = 2*colPtrM[3];
840 colPtrV[4] = 2*colPtrM[4];
841 colPtrV[5] = 2*colPtrM[5];
842 colPtrV[6] = 2*colPtrM[6];
843 colPtrV[7] = 2*colPtrM[7];
844 colPtrV[8] = 2*colPtrM[8];
846 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],16),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
848 colPtrV[0] = 2*colPtrM[0]+1;
849 colPtrV[1] = 2*colPtrM[1]+1;
850 colPtrV[2] = 2*colPtrM[2]+1;
851 colPtrV[3] = 2*colPtrM[3]+1;
852 colPtrV[4] = 2*colPtrM[4]+1;
853 colPtrV[5] = 2*colPtrM[5]+1;
854 colPtrV[6] = 2*colPtrM[6]+1;
855 colPtrV[7] = 2*colPtrM[7]+1;
856 colPtrV[8] = 2*colPtrM[8]+1;
858 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[1],17),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
861 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
862 valPtrM[0] = -0.046875;
863 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,1);
864 valPtrM[1] = 0.015625;
865 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,2);
866 valPtrM[2] = -0.046875;
867 colPtrM[3] = (*coarseElementMDOFs)(coarseElement,3);
868 valPtrM[3] = 0.140625;
869 colPtrM[4] = (*coarseElementMDOFs)(coarseElement,4);
870 valPtrM[4] = -0.09375;
871 colPtrM[5] = (*coarseElementMDOFs)(coarseElement,5);
872 valPtrM[5] = -0.09375;
873 colPtrM[6] = (*coarseElementMDOFs)(coarseElement,6);
874 valPtrM[6] = 0.28125;
875 colPtrM[7] = (*coarseElementMDOFs)(coarseElement,7);
876 valPtrM[7] = 0.28125;
877 colPtrM[8] = (*coarseElementMDOFs)(coarseElement,8);
881 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[2],8),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
883 colPtrV[0] = 2*colPtrM[0];
884 colPtrV[1] = 2*colPtrM[1];
885 colPtrV[2] = 2*colPtrM[2];
886 colPtrV[3] = 2*colPtrM[3];
887 colPtrV[4] = 2*colPtrM[4];
888 colPtrV[5] = 2*colPtrM[5];
889 colPtrV[6] = 2*colPtrM[6];
890 colPtrV[7] = 2*colPtrM[7];
891 colPtrV[8] = 2*colPtrM[8];
893 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],16),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
895 colPtrV[0] = 2*colPtrM[0]+1;
896 colPtrV[1] = 2*colPtrM[1]+1;
897 colPtrV[2] = 2*colPtrM[2]+1;
898 colPtrV[3] = 2*colPtrM[3]+1;
899 colPtrV[4] = 2*colPtrM[4]+1;
900 colPtrV[5] = 2*colPtrM[5]+1;
901 colPtrV[6] = 2*colPtrM[6]+1;
902 colPtrV[7] = 2*colPtrM[7]+1;
903 colPtrV[8] = 2*colPtrM[8]+1;
905 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[2],17),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
908 colPtrM[0] = (*coarseElementMDOFs)(coarseElement,0);
909 valPtrM[0] = 0.015625;
910 colPtrM[1] = (*coarseElementMDOFs)(coarseElement,1);
911 valPtrM[1] = -0.046875;
912 colPtrM[2] = (*coarseElementMDOFs)(coarseElement,2);
913 valPtrM[2] = 0.140625;
914 colPtrM[3] = (*coarseElementMDOFs)(coarseElement,3);
915 valPtrM[3] = -0.046875;
916 colPtrM[4] = (*coarseElementMDOFs)(coarseElement,4);
917 valPtrM[4] = -0.09375;
918 colPtrM[5] = (*coarseElementMDOFs)(coarseElement,5);
919 valPtrM[5] = 0.28125;
920 colPtrM[6] = (*coarseElementMDOFs)(coarseElement,6);
921 valPtrM[6] = 0.28125;
922 colPtrM[7] = (*coarseElementMDOFs)(coarseElement,7);
923 valPtrM[7] = -0.09375;
924 colPtrM[8] = (*coarseElementMDOFs)(coarseElement,8);
928 PM->insertGlobalValues((*fineElementMDOFs)(fineElement[3],8),colPtrM.
view(0,nnz),valPtrM.view(0,nnz));
930 colPtrV[0] = 2*colPtrM[0];
931 colPtrV[1] = 2*colPtrM[1];
932 colPtrV[2] = 2*colPtrM[2];
933 colPtrV[3] = 2*colPtrM[3];
934 colPtrV[4] = 2*colPtrM[4];
935 colPtrV[5] = 2*colPtrM[5];
936 colPtrV[6] = 2*colPtrM[6];
937 colPtrV[7] = 2*colPtrM[7];
938 colPtrV[8] = 2*colPtrM[8];
940 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],16),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
942 colPtrV[0] = 2*colPtrM[0]+1;
943 colPtrV[1] = 2*colPtrM[1]+1;
944 colPtrV[2] = 2*colPtrM[2]+1;
945 colPtrV[3] = 2*colPtrM[3]+1;
946 colPtrV[4] = 2*colPtrM[4]+1;
947 colPtrV[5] = 2*colPtrM[5]+1;
948 colPtrV[6] = 2*colPtrM[6]+1;
949 colPtrV[7] = 2*colPtrM[7]+1;
950 colPtrV[8] = 2*colPtrM[8]+1;
952 PV->insertGlobalValues((*fineElementVDOFs)(fineElement[3],17),colPtrV.
view(0,nnz),valPtrM.view(0,nnz));
956 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,0);
958 colPtrP[1] = (*coarseElementPDOFs)(coarseElement,1);
960 colPtrP[2] = (*coarseElementPDOFs)(coarseElement,2);
962 colPtrP[3] = (*coarseElementPDOFs)(coarseElement,3);
966 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[0],2),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
969 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,1);
971 colPtrP[1] = (*coarseElementPDOFs)(coarseElement,2);
975 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[1],2),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
978 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,2);
982 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[3],2),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
985 colPtrP[0] = (*coarseElementPDOFs)(coarseElement,2);
987 colPtrP[1] = (*coarseElementPDOFs)(coarseElement,3);
991 PP->insertGlobalValues((*fineElementPDOFs)(fineElement[3],3),colPtrP.
view(0,nnz),valPtrP.view(0,nnz));
995 if ((coarseElement+1) % (nCoarseElements) == 0)
997 fineElement[0] = fineElement[3]+1;
998 fineElement[1] = fineElement[0]+1;
999 fineElement[2] = fineElement[0]+nFineElements;
1000 fineElement[3] = fineElement[2]+1;
1004 fineElement[0] = fineElement[1]+1;
1005 fineElement[1] = fineElement[0]+1;
1006 fineElement[2] = fineElement[3]+1;
1007 fineElement[3] = fineElement[2]+1;
1014 for (GO VRow = 0; VRow < fNV; VRow++)
1019 PV->getGlobalRowView(VRow,colPtr,valPtr);
1022 P->insertGlobalValues(VRow,colPtr,valPtr);
1027 for (GO PRow = 0; PRow < fNP; PRow++)
1033 PP->getGlobalRowView(PRow,colPtr,valPtr);
1036 for (LO jj=0; jj<colPtr.
size(); jj++)
1038 newColPtr[jj] += colPtr[jj];
1042 P->insertGlobalValues(PRow+fNV,newColPtr.view(0,colPtr.
size()),valPtr);
1047 for (GO MRow = 0; MRow < fNM; MRow++)
1053 PM->getGlobalRowView(MRow,colPtr,valPtr);
1056 for (LO jj=0; jj<colPtr.
size(); jj++)
1058 newColPtr[jj] += colPtr[jj];
1062 P->insertGlobalValues(MRow+fNV+fNP,newColPtr.view(0,colPtr.
size()),valPtr);
1073 PV->fillComplete(colMapforPV,rowMapforPV);
1074 PP->fillComplete(colMapforPP,rowMapforPP);
1075 PM->fillComplete(colMapforPM,rowMapforPM);
1079 Set(coarseLevel,
"PV",PV);
1080 Set(coarseLevel,
"PP",PP);
1081 Set(coarseLevel,
"PM",PM);
1082 Set(coarseLevel,
"P" ,P );
1084 Set(coarseLevel,
"NV",nV);
1085 Set(coarseLevel,
"NP",nP);
1086 Set(coarseLevel,
"NM",nM);
1093 #define MUELU_GEOINTERPFACTORY_SHORT
1094 #endif // MUELU_GEOINTERPFACTORY_DEF_HPP
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.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
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.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
virtual ~GeoInterpFactory()
Destructor.
const RCP< const FactoryManagerBase > GetFactoryManager()
returns the current factory manager
Description of what is happening (more verbose)
GeoInterpFactory()
Constructor.
ArrayView< T > view(size_type lowerOffset, size_type size) const