46#ifndef MUELU_IPCFACTORY_DEF_HPP
47#define MUELU_IPCFACTORY_DEF_HPP
49#include <Xpetra_Matrix.hpp>
50#include <Xpetra_IO.hpp>
59#include "MueLu_PerfUtils.hpp"
60#include "MueLu_Utilities.hpp"
62#include "Teuchos_ScalarTraits.hpp"
70#include "Intrepid2_HGRAD_LINE_C1_FEM.hpp"
71#include "Intrepid2_HGRAD_LINE_Cn_FEM.hpp"
76#include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp"
78#include "Intrepid2_HGRAD_QUAD_Cn_FEM.hpp"
93#define MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(level,ename,entry) \
94 {if (level.IsRequested(ename,this) || level.GetKeepFlag(ename,this) != 0) this->Set(level,ename,entry);}
102namespace MueLuIntrepid {
103inline std::string
tolower(
const std::string & str) {
104 std::string data(str);
105 std::transform(data.begin(), data.end(), data.begin(),
::tolower);
111 template<
class Basis,
class LOFieldContainer,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
113 std::vector<std::vector<LocalOrdinal> > &seeds,
114 const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &rowMap,
115 const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &columnMap)
127 shards::CellTopology cellTopo = basis->getBaseCellTopology();
128 int spaceDim = cellTopo.getDimension();
130 seeds.resize(spaceDim + 1);
134 LocalOrdinal lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
135 GlobalOrdinal go_invalid = Teuchos::OrdinalTraits<GO>::invalid();
138 std::vector<std::set<LocalOrdinal>> seedSets(spaceDim+1);
140 int numCells = elementToNodeMap.extent(0);
141 auto elementToNodeMap_host = Kokkos::create_mirror_view(elementToNodeMap);
142 Kokkos::deep_copy(elementToNodeMap_host,elementToNodeMap);
143 for (
int cellOrdinal=0; cellOrdinal<numCells; cellOrdinal++)
145 for (
int d=0; d<=spaceDim; d++)
147 int subcellCount = cellTopo.getSubcellCount(d);
148 for (
int subcord=0; subcord<subcellCount; subcord++)
150 int dofCount = basis->getDofCount(d,subcord);
151 if (dofCount == 0)
continue;
153 GO leastGlobalDofOrdinal = go_invalid;
154 LO LID_leastGlobalDofOrdinal = lo_invalid;
155 for (
int basisOrdinalOrdinal=0; basisOrdinalOrdinal<dofCount; basisOrdinalOrdinal++)
157 int basisOrdinal = basis->getDofOrdinal(d,subcord,basisOrdinalOrdinal);
158 int colLID = elementToNodeMap_host(cellOrdinal,basisOrdinal);
159 if (colLID != Teuchos::OrdinalTraits<LO>::invalid())
163 if (rowLID != lo_invalid)
165 if ((leastGlobalDofOrdinal == go_invalid) || (colGID < leastGlobalDofOrdinal))
168 leastGlobalDofOrdinal = colGID;
169 LID_leastGlobalDofOrdinal = rowLID;
174 if (leastGlobalDofOrdinal != go_invalid)
176 seedSets[d].insert(LID_leastGlobalDofOrdinal);
181 for (
int d=0; d<=spaceDim;d++)
183 seeds[d] = std::vector<LocalOrdinal>(seedSets[d].begin(),seedSets[d].end());
194template<
class Scalar,
class KokkosExecutionSpace>
195Teuchos::RCP< Intrepid2::Basis<KokkosExecutionSpace,Scalar,Scalar> >
BasisFactory(
const std::string & name,
int & degree)
199 string myerror(
"IntrepidBasisFactory: cannot parse string name '"+name+
"'");
204 size_t pos1 = name.find_first_of(
" _");
205 if(pos1==0)
throw std::runtime_error(myerror);
206 string deriv =
tolower(name.substr(0,pos1));
207 if(deriv!=
"hgrad" && deriv!=
"hcurl" && deriv!=
"hdiv")
throw std::runtime_error(myerror);
211 size_t pos2 = name.find_first_of(
" _",pos1);
212 if(pos2==0)
throw std::runtime_error(myerror);
213 string el =
tolower(name.substr(pos1,pos2-pos1));
214 if(el!=
"hex" && el!=
"line" && el!=
"poly" && el!=
"pyr" && el!=
"quad" && el!=
"tet" && el!=
"tri" && el!=
"wedge")
throw std::runtime_error(myerror);
218 string poly =
tolower(name.substr(pos2,1));
219 if(poly!=
"c" && poly!=
"i")
throw std::runtime_error(myerror);
223 degree=std::stoi(name.substr(pos2,1));
224 if(degree<=0)
throw std::runtime_error(myerror);
227 if(deriv==
"hgrad" && el==
"quad" && poly==
"c"){
228 if(degree==1)
return rcp(
new Intrepid2::Basis_HGRAD_QUAD_C1_FEM<KokkosExecutionSpace,Scalar,Scalar>());
229 else return rcp(
new Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar>(degree,Intrepid2::POINTTYPE_EQUISPACED));
231 else if(deriv==
"hgrad" && el==
"line" && poly==
"c"){
232 if(degree==1)
return rcp(
new Intrepid2::Basis_HGRAD_LINE_C1_FEM<KokkosExecutionSpace,Scalar,Scalar>());
233 else return rcp(
new Intrepid2::Basis_HGRAD_LINE_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar>(degree,Intrepid2::POINTTYPE_EQUISPACED));
237 throw std::runtime_error(myerror);
238 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
248template<
class Scalar,
class KokkosDeviceType>
249void IntrepidGetP1NodeInHi(
const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space,Scalar,Scalar> >&hi_basis,
250 std::vector<size_t> & lo_node_in_hi,
251 Kokkos::DynRankView<Scalar,KokkosDeviceType> & hi_DofCoords) {
253 typedef typename KokkosDeviceType::execution_space KokkosExecutionSpace;
255 size_t degree = hi_basis->getDegree();
256 lo_node_in_hi.resize(0);
258 if(!rcp_dynamic_cast<Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar> >(hi_basis).is_null()) {
260 lo_node_in_hi.insert(lo_node_in_hi.end(),{0,degree, (degree+1)*(degree+1)-1, degree*(degree+1)});
262 else if(!rcp_dynamic_cast<Intrepid2::Basis_HGRAD_LINE_Cn_FEM<KokkosExecutionSpace,Scalar,Scalar> >(hi_basis).is_null()) {
264 lo_node_in_hi.insert(lo_node_in_hi.end(),{0,degree});
267 throw std::runtime_error(
"IntrepidPCoarsenFactory: Unknown element type");
270 Kokkos::resize(hi_DofCoords,hi_basis->getCardinality(),hi_basis->getBaseCellTopology().getDimension());
271 hi_basis->getDofCoords(hi_DofCoords);
283template<
class LocalOrdinal,
class GlobalOrdinal,
class Node,
class LOFieldContainer>
285 RCP<
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & hi_columnMap,
286 LOFieldContainer & lo_elemToHiRepresentativeNode) {
292 size_t numElem = hi_elemToNode.extent(0);
293 size_t lo_nperel = candidates.size();
294 Kokkos::resize(lo_elemToHiRepresentativeNode,numElem, lo_nperel);
296 auto lo_elemToHiRepresentativeNode_host = Kokkos::create_mirror_view(lo_elemToHiRepresentativeNode);
297 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
298 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
299 for(
size_t i=0; i<numElem; i++)
300 for(
size_t j=0; j<lo_nperel; j++) {
301 if(candidates[j].size() == 1)
302 lo_elemToHiRepresentativeNode_host(i,j) = hi_elemToNode_host(i,candidates[j][0]);
305 std::vector<GO> GID(candidates[j].size());
306 for(
size_t k=0; k<(size_t)candidates[j].size(); k++)
307 GID[k] = hi_columnMap->getGlobalElement(hi_elemToNode_host(i,candidates[j][k]));
310 size_t which = std::distance(GID.begin(),std::min_element(GID.begin(),GID.end()));
313 lo_elemToHiRepresentativeNode_host(i,j) = hi_elemToNode_host(i,candidates[j][which]);
316 Kokkos::deep_copy(lo_elemToHiRepresentativeNode, lo_elemToHiRepresentativeNode_host);
329template <
class LocalOrdinal,
class LOFieldContainer>
331 const std::vector<bool> & hi_nodeIsOwned,
332 const LOFieldContainer & lo_elemToHiRepresentativeNode,
333 LOFieldContainer & lo_elemToNode,
334 std::vector<bool> & lo_nodeIsOwned,
335 std::vector<LocalOrdinal> & hi_to_lo_map,
336 int & lo_numOwnedNodes) {
340 size_t numElem = hi_elemToNode.extent(0);
341 size_t hi_numNodes = hi_nodeIsOwned.size();
342 size_t lo_nperel = lo_elemToHiRepresentativeNode.extent(1);
343 Kokkos::resize(lo_elemToNode,numElem, lo_nperel);
346 auto lo_elemToHiRepresentativeNode_host = Kokkos::create_mirror_view(lo_elemToHiRepresentativeNode);
347 Kokkos::deep_copy(lo_elemToHiRepresentativeNode_host, lo_elemToHiRepresentativeNode);
348 std::vector<bool> is_low_order(hi_numNodes,
false);
349 for(
size_t i=0; i<numElem; i++)
350 for(
size_t j=0; j<lo_nperel; j++) {
351 LO
id = lo_elemToHiRepresentativeNode_host(i,j);
352 is_low_order[id] =
true;
357 size_t lo_numNodes=0;
358 hi_to_lo_map.resize(hi_numNodes,Teuchos::OrdinalTraits<LO>::invalid());
360 for(
size_t i=0; i<hi_numNodes; i++)
361 if(is_low_order[i]) {
362 hi_to_lo_map[i] = lo_numNodes;
364 if(hi_nodeIsOwned[i]) lo_numOwnedNodes++;
368 lo_nodeIsOwned.resize(lo_numNodes,
false);
369 for(
size_t i=0; i<hi_numNodes; i++) {
370 if(is_low_order[i] && hi_nodeIsOwned[i])
371 lo_nodeIsOwned[hi_to_lo_map[i]]=
true;
375 auto lo_elemToNode_host = Kokkos::create_mirror_view(lo_elemToNode);
376 for(
size_t i=0; i<numElem; i++)
377 for(
size_t j=0; j<lo_nperel; j++)
378 lo_elemToNode_host(i,j) = hi_to_lo_map[lo_elemToHiRepresentativeNode_host(i,j)];
383 bool map_ordering_test_passed=
true;
384 for(
size_t i=0; i<lo_numNodes-1; i++)
385 if(!lo_nodeIsOwned[i] && lo_nodeIsOwned[i+1])
386 map_ordering_test_passed=
false;
388 if(!map_ordering_test_passed)
389 throw std::runtime_error(
"MueLu::MueLuIntrepid::BuildLoElemToNodeViaRepresentatives failed map ordering test");
390 Kokkos::deep_copy(lo_elemToNode, lo_elemToNode_host);
406template <
class LocalOrdinal,
class LOFieldContainer>
408 const std::vector<bool> & hi_nodeIsOwned,
409 const std::vector<size_t> & lo_node_in_hi,
410 const Teuchos::ArrayRCP<const int> & hi_isDirichlet,
411 LOFieldContainer & lo_elemToNode,
412 std::vector<bool> & lo_nodeIsOwned,
413 std::vector<LocalOrdinal> & hi_to_lo_map,
414 int & lo_numOwnedNodes) {
417 LocalOrdinal LOINVALID = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
420 size_t numElem = hi_elemToNode.extent(0);
421 size_t hi_numNodes = hi_nodeIsOwned.size();
423 size_t lo_nperel = lo_node_in_hi.size();
424 Kokkos::resize(lo_elemToNode,numElem, lo_nperel);
427 std::vector<bool> is_low_order(hi_numNodes,
false);
428 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
429 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
430 auto lo_elemToNode_host = Kokkos::create_mirror_view(lo_elemToNode);
431 for(
size_t i=0; i<numElem; i++)
432 for(
size_t j=0; j<lo_nperel; j++) {
433 LO lid = hi_elemToNode_host(i,lo_node_in_hi[j]);
436 if(hi_isDirichlet[lid])
437 lo_elemToNode_host(i,j) = LOINVALID;
439 lo_elemToNode_host(i,j) = lid;
440 is_low_order[hi_elemToNode_host(i,lo_node_in_hi[j])] =
true;
446 size_t lo_numNodes=0;
447 hi_to_lo_map.resize(hi_numNodes,Teuchos::OrdinalTraits<LO>::invalid());
449 for(
size_t i=0; i<hi_numNodes; i++)
450 if(is_low_order[i]) {
451 hi_to_lo_map[i] = lo_numNodes;
453 if(hi_nodeIsOwned[i]) lo_numOwnedNodes++;
457 lo_nodeIsOwned.resize(lo_numNodes,
false);
458 for(
size_t i=0; i<hi_numNodes; i++) {
459 if(is_low_order[i] && hi_nodeIsOwned[i])
460 lo_nodeIsOwned[hi_to_lo_map[i]]=
true;
464 for(
size_t i=0; i<numElem; i++)
465 for(
size_t j=0; j<lo_nperel; j++) {
466 if(lo_elemToNode_host(i,j) != LOINVALID)
467 lo_elemToNode_host(i,j) = hi_to_lo_map[lo_elemToNode_host(i,j)];
469 Kokkos::deep_copy(lo_elemToNode, lo_elemToNode_host);
473 bool map_ordering_test_passed=
true;
474 for(
size_t i=0; i<lo_numNodes-1; i++)
475 if(!lo_nodeIsOwned[i] && lo_nodeIsOwned[i+1])
476 map_ordering_test_passed=
false;
478 if(!map_ordering_test_passed)
479 throw std::runtime_error(
"MueLu::MueLuIntrepid::BuildLoElemToNode failed map ordering test");
492 template <
class LocalOrdinal,
class GlobalOrdinal,
class Node>
493 void GenerateColMapFromImport(
const Xpetra::Import<LocalOrdinal,GlobalOrdinal, Node> & hi_importer,
const std::vector<LocalOrdinal> &hi_to_lo_map,
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> & lo_domainMap,
const size_t & lo_columnMapLength, RCP<
const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & lo_columnMap) {
497 typedef Xpetra::Map<LO,GO,NO> Map;
498 typedef Xpetra::Vector<GO,LO,GO,NO> GOVector;
500 GO go_invalid = Teuchos::OrdinalTraits<GO>::invalid();
501 LO lo_invalid = Teuchos::OrdinalTraits<LO>::invalid();
503 RCP<const Map> hi_domainMap = hi_importer.getSourceMap();
504 RCP<const Map> hi_columnMap = hi_importer.getTargetMap();
510 RCP<GOVector> dvec = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(hi_domainMap);
512 ArrayRCP<GO> dvec_data = dvec->getDataNonConst(0);
513 for(
size_t i=0; i<hi_domainMap->getLocalNumElements(); i++) {
514 if(hi_to_lo_map[i]!=lo_invalid) dvec_data[i] = lo_domainMap.getGlobalElement(hi_to_lo_map[i]);
515 else dvec_data[i] = go_invalid;
520 RCP<GOVector> cvec = Xpetra::VectorFactory<GO, LO, GO, NO>::Build(hi_columnMap,
true);
521 cvec->doImport(*dvec,hi_importer,Xpetra::ADD);
525 Array<GO> lo_col_data(lo_columnMapLength);
527 ArrayRCP<GO> cvec_data = cvec->getDataNonConst(0);
528 for(
size_t i=0,idx=0; i<hi_columnMap->getLocalNumElements(); i++) {
529 if(hi_to_lo_map[i]!=lo_invalid) {
530 lo_col_data[idx] = cvec_data[i];
536 lo_columnMap = Xpetra::MapFactory<LO,GO,NO>::Build(lo_domainMap.lib(),Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),lo_col_data(),lo_domainMap.getIndexBase(),lo_domainMap.getComm());
547template<
class Basis,
class SCFieldContainer>
548void GenerateRepresentativeBasisNodes(
const Basis & basis,
const SCFieldContainer & ReferenceNodeLocations,
const double threshold,std::vector<std::vector<size_t> > & representative_node_candidates) {
549 typedef SCFieldContainer FC;
550 typedef typename FC::data_type SC;
553 size_t numFieldsHi = ReferenceNodeLocations.extent(0);
555 size_t numFieldsLo = basis.getCardinality();
557 FC LoValues(
"LoValues",numFieldsLo,numFieldsHi);
559 basis.getValues(LoValues, ReferenceNodeLocations , Intrepid2::OPERATOR_VALUE);
564 printf(
"** LoValues[%d,%d] **\n",(
int)numFieldsLo,(
int)numFieldsHi);
565 for(
size_t i=0; i<numFieldsLo; i++) {
566 for(
size_t j=0; j<numFieldsHi; j++)
567 printf(
"%6.4e ",LoValues(i,j));
570 printf(
"**************\n");fflush(stdout);
573 representative_node_candidates.resize(numFieldsLo);
574 auto LoValues_host = Kokkos::create_mirror_view(LoValues);
575 Kokkos::deep_copy(LoValues_host, LoValues);
576 for(
size_t i=0; i<numFieldsLo; i++) {
578 typename Teuchos::ScalarTraits<SC>::magnitudeType vmax = Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<SC>::magnitudeType>::zero();
579 for(
size_t j=0; j<numFieldsHi; j++)
580 vmax = std::max(vmax,Teuchos::ScalarTraits<SC>::magnitude(LoValues_host(i,j)));
583 for(
size_t j=0; j<numFieldsHi; j++) {
584 if(Teuchos::ScalarTraits<SC>::magnitude(vmax - LoValues_host(i,j)) < threshold*vmax)
585 representative_node_candidates[i].push_back(j);
590 for(
size_t i=0; i<numFieldsLo; i++)
591 if(!representative_node_candidates[i].size())
592 throw std::runtime_error(
"ERROR: GenerateRepresentativeBasisNodes: No candidates found!");
605template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
607 const std::vector<bool> & hi_nodeIsOwned,
609 const std::vector<size_t> &lo_node_in_hi,
610 const Basis &lo_basis,
611 const std::vector<LocalOrdinal> & hi_to_lo_map,
612 const Teuchos::RCP<const Map> & lo_colMap,
613 const Teuchos::RCP<const Map> & lo_domainMap,
614 const Teuchos::RCP<const Map> & hi_map,
615 Teuchos::RCP<Matrix>& P)
const{
618 size_t numFieldsHi = hi_elemToNode.extent(1);
619 size_t numFieldsLo = lo_basis.getCardinality();
620 LocalOrdinal LOINVALID = Teuchos::OrdinalTraits<LocalOrdinal>::invalid();
621 FC LoValues_at_HiDofs(
"LoValues_at_HiDofs",numFieldsLo,numFieldsHi);
622 lo_basis.getValues(LoValues_at_HiDofs, hi_DofCoords, Intrepid2::OPERATOR_VALUE);
623 auto LoValues_at_HiDofs_host = Kokkos::create_mirror_view(LoValues_at_HiDofs);
624 Kokkos::deep_copy(LoValues_at_HiDofs_host, LoValues_at_HiDofs);
627 typedef typename Teuchos::ScalarTraits<SC>::halfPrecision SClo;
628 typedef typename Teuchos::ScalarTraits<SClo>::magnitudeType MT;
629 MT effective_zero = Teuchos::ScalarTraits<MT>::eps();
632 P = rcp(
new CrsMatrixWrap(hi_map,lo_colMap,numFieldsHi));
633 RCP<CrsMatrix> Pcrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
636 size_t Nelem=hi_elemToNode.extent(0);
637 std::vector<bool> touched(hi_map->getLocalNumElements(),
false);
638 Teuchos::Array<GO> col_gid(1);
639 Teuchos::Array<SC> val(1);
640 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
641 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
642 for(
size_t i=0; i<Nelem; i++) {
643 for(
size_t j=0; j<numFieldsHi; j++) {
644 LO row_lid = hi_elemToNode_host(i,j);
645 GO row_gid = hi_map->getGlobalElement(row_lid);
646 if(hi_nodeIsOwned[row_lid] && !touched[row_lid]) {
647 for(
size_t k=0; k<numFieldsLo; k++) {
649 LO col_lid = hi_to_lo_map[hi_elemToNode_host(i,lo_node_in_hi[k])];
650 if(col_lid==LOINVALID)
continue;
652 col_gid[0] = {lo_colMap->getGlobalElement(col_lid)};
653 val[0] = LoValues_at_HiDofs_host(k,j);
656 if(Teuchos::ScalarTraits<SC>::magnitude(val[0]) >= effective_zero)
657 P->insertGlobalValues(row_gid,col_gid(),val());
659 touched[row_lid]=
true;
663 P->fillComplete(lo_domainMap,hi_map);
667template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
669 const std::vector<bool> & hi_nodeIsOwned,
672 const Basis &lo_basis,
673 const std::vector<LocalOrdinal> & hi_to_lo_map,
674 const Teuchos::RCP<const Map> & lo_colMap,
675 const Teuchos::RCP<const Map> & lo_domainMap,
676 const Teuchos::RCP<const Map> & hi_map,
677 Teuchos::RCP<Matrix>& P)
const{
680 size_t numFieldsHi = hi_elemToNode.extent(1);
681 size_t numFieldsLo = lo_basis.getCardinality();
682 FC LoValues_at_HiDofs(
"LoValues_at_HiDofs",numFieldsLo,numFieldsHi);
683 lo_basis.getValues(LoValues_at_HiDofs, hi_DofCoords, Intrepid2::OPERATOR_VALUE);
684 auto LoValues_at_HiDofs_host = Kokkos::create_mirror_view(LoValues_at_HiDofs);
685 auto hi_elemToNode_host = Kokkos::create_mirror_view(hi_elemToNode);
686 auto lo_elemToHiRepresentativeNode_host = Kokkos::create_mirror_view(lo_elemToHiRepresentativeNode);
687 Kokkos::deep_copy(LoValues_at_HiDofs_host, LoValues_at_HiDofs);
688 Kokkos::deep_copy(hi_elemToNode_host, hi_elemToNode);
689 Kokkos::deep_copy(lo_elemToHiRepresentativeNode_host, lo_elemToHiRepresentativeNode);
692 typedef typename Teuchos::ScalarTraits<SC>::halfPrecision SClo;
693 typedef typename Teuchos::ScalarTraits<SClo>::magnitudeType MT;
694 MT effective_zero = Teuchos::ScalarTraits<MT>::eps();
697 P = rcp(
new CrsMatrixWrap(hi_map,lo_colMap,numFieldsHi));
698 RCP<CrsMatrix> Pcrs = rcp_dynamic_cast<CrsMatrixWrap>(P)->getCrsMatrix();
701 size_t Nelem=hi_elemToNode.extent(0);
702 std::vector<bool> touched(hi_map->getLocalNumElements(),
false);
703 Teuchos::Array<GO> col_gid(1);
704 Teuchos::Array<SC> val(1);
705 for(
size_t i=0; i<Nelem; i++) {
706 for(
size_t j=0; j<numFieldsHi; j++) {
707 LO row_lid = hi_elemToNode_host(i,j);
708 GO row_gid = hi_map->getGlobalElement(row_lid);
709 if(hi_nodeIsOwned[row_lid] && !touched[row_lid]) {
710 for(
size_t k=0; k<numFieldsLo; k++) {
712 LO col_lid = hi_to_lo_map[lo_elemToHiRepresentativeNode_host(i,k)];
713 col_gid[0] = {lo_colMap->getGlobalElement(col_lid)};
714 val[0] = LoValues_at_HiDofs_host(k,j);
717 if(Teuchos::ScalarTraits<SC>::magnitude(val[0]) >= effective_zero)
718 P->insertGlobalValues(row_gid,col_gid(),val());
720 touched[row_lid]=
true;
724 P->fillComplete(lo_domainMap,hi_map);
728 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
730 RCP<ParameterList> validParamList = rcp(
new ParameterList());
732#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
735#undef SET_VALID_ENTRY
737 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A used during the prolongator smoothing process");
739 validParamList->set< RCP<const FactoryBase> >(
"Nullspace", Teuchos::null,
"Generating factory of the nullspace");
740 validParamList->set< RCP<const FactoryBase> >(
"pcoarsen: element to node map", Teuchos::null,
"Generating factory of the element to node map");
741 return validParamList;
745 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
747 Input(fineLevel,
"A");
748 Input(fineLevel,
"pcoarsen: element to node map");
749 Input(fineLevel,
"Nullspace");
753 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
755 return BuildP(fineLevel, coarseLevel);
759 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
763 const std::string prefix =
"MueLu::IntrepidPCoarsenFactory(" + levelIDs +
"): ";
766 typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> FCi;
767 typedef Kokkos::DynRankView<double,typename Node::device_type> FC;
772 Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>& Acrs =
dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>&
>(*A);
781 std::vector<LocalOrdinal> A_dirichletRows;
788 RCP<ParameterList> APparams = rcp(
new ParameterList);
789 if (coarseLevel.
IsAvailable(
"AP reuse data",
this)) {
792 APparams = coarseLevel.
Get< RCP<ParameterList> >(
"AP reuse data",
this);
794 if (APparams->isParameter(
"graph"))
795 finalP = APparams->get< RCP<Matrix> >(
"graph");
804 int lo_degree, hi_degree;
809 GetOStream(
Statistics1) <<
"P-Coarsening from basis "<<pL.get<std::string>(
"pcoarsen: hi basis")<<
" to "<<pL.get<std::string>(
"pcoarsen: lo basis") <<std::endl;
813 const Teuchos::RCP<FCi> Pn_elemToNode =
Get<Teuchos::RCP<FCi> >(fineLevel,
"pcoarsen: element to node map");
820 RCP<const Map> rowMap = A->getRowMap();
821 RCP<const Map> colMap = Acrs.getColMap();
822 RCP<const Map> domainMap = A->getDomainMap();
823 int NumProc = rowMap->getComm()->getSize();
824 assert(rowMap->isSameAs(*domainMap));
825 std::vector<bool> Pn_nodeIsOwned(colMap->getLocalNumElements(),
false);
826 LO num_owned_rows = 0;
827 for(
size_t i=0; i<rowMap->getLocalNumElements(); i++) {
828 if(rowMap->getGlobalElement(i) == colMap->getGlobalElement(i)) {
829 Pn_nodeIsOwned[i] =
true;
836 Teuchos::RCP<FCi> P1_elemToNode = rcp(
new FCi());
838 std::vector<bool> P1_nodeIsOwned;
839 int P1_numOwnedNodes;
840 std::vector<LO> hi_to_lo_map;
843 std::vector<size_t> lo_node_in_hi;
846 FCi lo_elemToHiRepresentativeNode;
849 RCP<Xpetra::Vector<int,LocalOrdinal,GlobalOrdinal,Node> > hi_isDirichletRow, hi_isDirichletCol;
853 printf(
"[%d] isDirichletRow = ",A->getRowMap()->getComm()->getRank());
854 for(
size_t i=0;i<hi_isDirichletRow->getMap()->getLocalNumElements(); i++)
855 printf(
"%d ",hi_isDirichletRow->getData(0)[i]);
857 printf(
"[%d] isDirichletCol = ",A->getRowMap()->getComm()->getRank());
858 for(
size_t i=0;i<hi_isDirichletCol->getMap()->getLocalNumElements(); i++)
859 printf(
"%d ",hi_isDirichletCol->getData(0)[i]);
870 MueLuIntrepid::BuildLoElemToNode(*Pn_elemToNode,Pn_nodeIsOwned,lo_node_in_hi,hi_isDirichletCol->getData(0),*P1_elemToNode,P1_nodeIsOwned,hi_to_lo_map,P1_numOwnedNodes);
871 assert(hi_to_lo_map.size() == colMap->getLocalNumElements());
875 double threshold = 1e-10;
876 std::vector<std::vector<size_t> > candidates;
877 Kokkos::resize(hi_DofCoords,hi_basis->getCardinality(),hi_basis->getBaseCellTopology().getDimension());
878 hi_basis->getDofCoords(hi_DofCoords);
891 RCP<const Map> P1_domainMap = MapFactory::Build(rowMap->lib(),Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(),P1_numOwnedNodes,rowMap->getIndexBase(),rowMap->getComm());
895 RCP<const Map> P1_colMap;
896 if(NumProc==1) P1_colMap = P1_domainMap;
912 RCP<MultiVector> coarseNullspace = MultiVectorFactory::Build(P1_domainMap, fineNullspace->getNumVectors());
913 finalP->apply(*fineNullspace,*coarseNullspace,Teuchos::TRANS);
914 Set(coarseLevel,
"Nullspace", coarseNullspace);
919 Set(coarseLevel,
"P", finalP);
921 APparams->set(
"graph", finalP);
925 RCP<ParameterList> params = rcp(
new ParameterList());
926 params->set(
"printLoadBalancingInfo",
true);
927 params->set(
"printCommInfo",
true);
938 Set(coarseLevel,
"R", R);
941 RCP<ParameterList> params = rcp(
new ParameterList());
942 params->set(
"printLoadBalancingInfo",
true);
943 params->set(
"printCommInfo",
true);
#define SET_VALID_ENTRY(name)
#define MUELU_LEVEL_SET_IF_REQUESTED_OR_KEPT(level, ename, entry)
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Kokkos::DynRankView< LocalOrdinal, typename Node::device_type > LOFieldContainer
void GenerateLinearCoarsening_pn_kirby_to_p1(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const std::vector< size_t > &lo_node_in_hi, const Basis &lo_Basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
Kokkos::DynRankView< double, typename Node::device_type > SCFieldContainer
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Intrepid2::Basis< typename Node::device_type::execution_space, double, double > Basis
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void GenerateLinearCoarsening_pn_kirby_to_pm(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const LOFieldContainer &lo_elemToHiRepresentativeNode, const Basis &lo_basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
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.
int GetLevelID() const
Return level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
virtual const Teuchos::ParameterList & GetParameterList() const
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
static void FindDirichletRowsAndPropagateToCols(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, Teuchos::RCP< Xpetra::Vector< int, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &isDirichletRow, Teuchos::RCP< Xpetra::Vector< int, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &isDirichletCol)
static void ZeroDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, const std::vector< DefaultLocalOrdinal > &dirichletRows, Scalar replaceWith=Teuchos::ScalarTraits< Scalar >::zero())
static void FindDirichletRows(Teuchos::RCP< Xpetra::Matrix< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > > &A, std::vector< DefaultLocalOrdinal > &dirichletRows, bool count_twos_as_dirichlet=false)
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)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
void IntrepidGetP1NodeInHi(const Teuchos::RCP< Intrepid2::Basis< typename KokkosDeviceType::execution_space, Scalar, Scalar > > &hi_basis, std::vector< size_t > &lo_node_in_hi, Kokkos::DynRankView< Scalar, KokkosDeviceType > &hi_DofCoords)
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer &ReferenceNodeLocations, const double threshold, std::vector< std::vector< size_t > > &representative_node_candidates)
void GenerateColMapFromImport(const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &hi_importer, const std::vector< LocalOrdinal > &hi_to_lo_map, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &lo_domainMap, const size_t &lo_columnMapLength, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &lo_columnMap)
void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const LOFieldContainer &lo_elemToHiRepresentativeNode, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
std::string tolower(const std::string &str)
void BuildLoElemToNode(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const std::vector< size_t > &lo_node_in_hi, const Teuchos::ArrayRCP< const int > &hi_isDirichlet, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
void GenerateLoNodeInHiViaGIDs(const std::vector< std::vector< size_t > > &candidates, const LOFieldContainer &hi_elemToNode, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &hi_columnMap, LOFieldContainer &lo_elemToHiRepresentativeNode)
Teuchos::RCP< Intrepid2::Basis< KokkosExecutionSpace, Scalar, Scalar > > BasisFactory(const std::string &name, int °ree)
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Statistics1
Print more statistics.
@ Runtime0
One-liner description of what is happening.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.