46 #ifndef MUELU_TENTATIVEPFACTORY_KOKKOS_DEF_HPP
47 #define MUELU_TENTATIVEPFACTORY_KOKKOS_DEF_HPP
49 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
55 #include "MueLu_Aggregates_kokkos.hpp"
56 #include "MueLu_AmalgamationFactory_kokkos.hpp"
57 #include "MueLu_AmalgamationInfo_kokkos.hpp"
58 #include "MueLu_CoarseMapFactory_kokkos.hpp"
60 #include "MueLu_NullspaceFactory_kokkos.hpp"
61 #include "MueLu_PerfUtils.hpp"
63 #include "MueLu_Utilities_kokkos.hpp"
69 template<
class LocalOrdinal,
class View>
70 class ReduceMaxFunctor{
72 ReduceMaxFunctor(View view) : view_(view) { }
74 KOKKOS_INLINE_FUNCTION
80 KOKKOS_INLINE_FUNCTION
87 KOKKOS_INLINE_FUNCTION
96 template<
class LOType,
class GOType,
class SCType,
class DeviceType,
class NspType,
class aggRowsType,
class maxAggDofSizeType,
class agg2RowMapLOType,
class statusType,
class rowsType,
class rowsAuxType,
class colsAuxType,
class valsAuxType>
97 class LocalQRDecompFunctor {
103 typedef typename DeviceType::execution_space execution_space;
104 typedef Kokkos::ArithTraits<SC> ATS;
105 typedef typename ATS::magnitudeType Magnitude;
115 maxAggDofSizeType maxAggDofSize;
116 agg2RowMapLOType agg2RowMapLO;
117 statusType statusAtomic;
124 LocalQRDecompFunctor(NspType fineNS_, NspType coarseNS_, aggRowsType aggRows_, maxAggDofSizeType maxAggDofSize_, agg2RowMapLOType agg2RowMapLO_, statusType statusAtomic_, rowsType rows_, rowsAuxType rowsAux_, colsAuxType colsAux_, valsAuxType valsAux_,
bool doQRStep_) :
128 maxAggDofSize(maxAggDofSize_),
129 agg2RowMapLO(agg2RowMapLO_),
130 statusAtomic(statusAtomic_),
138 KOKKOS_INLINE_FUNCTION
140 auto agg = thread.league_rank();
143 auto aggSize = aggRows(agg+1) - aggRows(agg);
145 const SC one = ATS::one();
146 const SC two = one + one;
147 const SC zero = ATS::zero();
148 const auto zeroM = ATS::magnitude(zero);
151 int n = fineNS.extent(1);
159 shared_matrix r(thread.team_shmem(), m, n);
160 for (
int j = 0; j < n; j++)
161 for (
int k = 0; k < m; k++)
162 r(k,j) = fineNS(agg2RowMapLO(aggRows(agg)+k),j);
165 for (
int i = 0; i < m; i++) {
166 for (
int j = 0; j < n; j++)
167 printf(
" %5.3lf ", r(i,j));
173 shared_matrix q(thread.team_shmem(), m, m);
175 bool isSingular =
false;
179 for (
int i = 0; i < m; i++) {
180 for (
int j = 0; j < m; j++)
185 for (
int k = 0; k < n; k++) {
187 Magnitude s = zeroM, norm, norm_x;
188 for (
int i = k+1; i < m; i++)
189 s +=
pow(ATS::magnitude(r(i,k)), 2);
190 norm =
sqrt(
pow(ATS::magnitude(r(k,k)), 2) + s);
199 norm_x =
sqrt(
pow(ATS::magnitude(r(k,k)), 2) + s);
200 if (norm_x == zeroM) {
208 for (
int i = k; i < m; i++)
212 for (
int j = k+1; j < n; j++) {
215 for (
int i = k; i < m; i++)
216 si += r(i,k) * r(i,j);
217 for (
int i = k; i < m; i++)
218 r(i,j) -= two*si * r(i,k);
222 for (
int j = k; j < m; j++) {
225 for (
int i = k; i < m; i++)
226 si += r(i,k) * qt(i,j);
227 for (
int i = k; i < m; i++)
228 qt(i,j) -= two*si * r(i,k);
233 for (
int i = k+1; i < m; i++)
239 for (
int i = 0; i < m; i++)
240 for (
int j = 0; j < i; j++) {
248 for (
int j = 0; j < n; j++)
249 for (
int k = 0; k <= j; k++)
250 coarseNS(offset+k,j) = r(k,j);
253 statusAtomic(1) =
true;
288 for (
int j = 0; j < n; j++)
289 for (
int k = 0; k < n; k++)
291 coarseNS(offset+k,j) = r(k,j);
293 coarseNS(offset+k,j) = (k == j ? one : zero);
296 for (
int i = 0; i < m; i++)
297 for (
int j = 0; j < n; j++)
298 q(i,j) = (j == i ? one : zero);
302 for (
int j = 0; j < m; j++) {
303 LO localRow = agg2RowMapLO(aggRows(agg)+j);
304 size_t rowStart = rowsAux(localRow);
306 for (
int k = 0; k < n; k++) {
308 if (q(j,k) != zero) {
309 colsAux(rowStart+lnnz) = offset + k;
310 valsAux(rowStart+lnnz) = q(j,k);
314 rows(localRow+1) = lnnz;
320 for (
int i = 0; i < m; i++) {
321 for (
int j = 0; j < n; j++)
322 printf(
" %5.3lf ", coarseNS(i,j));
327 for (
int i = 0; i < aggSize; i++) {
328 for (
int j = 0; j < aggSize; j++)
329 printf(
" %5.3lf ", q(i,j));
343 for (
int j = 0; j < m; j++) {
344 LO localRow = agg2RowMapLO(aggRows(agg)+j);
345 size_t rowStart = rowsAux(localRow);
347 for (
int k = 0; k < n; k++) {
348 const SC qr_jk = fineNS(localRow,k);
351 colsAux(rowStart+lnnz) = offset + k;
352 valsAux(rowStart+lnnz) = qr_jk;
356 rows(localRow+1) = lnnz;
360 for (
int j = 0; j < n; j++)
361 coarseNS(offset+j,j) = one;
368 size_t team_shmem_size(
int )
const {
370 int m = maxAggDofSize;
371 int n = fineNS.extent(1);
372 return shared_matrix::shmem_size(m, n) +
373 shared_matrix::shmem_size(m, m);
381 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
382 RCP<const ParameterList> TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::GetValidParameterList()
const {
383 RCP<ParameterList> validParamList =
rcp(
new ParameterList());
385 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
388 #undef SET_VALID_ENTRY
390 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
391 validParamList->set< RCP<const FactoryBase> >(
"Aggregates", Teuchos::null,
"Generating factory of the aggregates");
392 validParamList->set< RCP<const FactoryBase> >(
"Nullspace", Teuchos::null,
"Generating factory of the nullspace");
393 validParamList->set< RCP<const FactoryBase> >(
"Scaled Nullspace", Teuchos::null,
"Generating factory of the scaled nullspace");
394 validParamList->set< RCP<const FactoryBase> >(
"UnAmalgamationInfo", Teuchos::null,
"Generating factory of UnAmalgamationInfo");
395 validParamList->set< RCP<const FactoryBase> >(
"CoarseMap", Teuchos::null,
"Generating factory of the coarse map");
396 validParamList->set< RCP<const FactoryBase> >(
"Coordinates", Teuchos::null,
"Generating factory of the coordinates");
399 ParameterList norecurse;
400 norecurse.disableRecursiveValidation();
401 validParamList->set<ParameterList> (
"matrixmatrix: kernel params", norecurse,
"MatrixMatrix kernel parameters");
403 return validParamList;
406 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
407 void TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::DeclareInput(Level& fineLevel, Level& )
const {
409 const ParameterList& pL = GetParameterList();
411 std::string nspName =
"Nullspace";
412 if(pL.isParameter(
"Nullspace name")) nspName = pL.get<std::string>(
"Nullspace name");
414 Input(fineLevel,
"A");
415 Input(fineLevel,
"Aggregates");
416 Input(fineLevel, nspName);
417 Input(fineLevel,
"UnAmalgamationInfo");
418 Input(fineLevel,
"CoarseMap");
419 if( fineLevel.GetLevelID() == 0 &&
421 pL.get<
bool>(
"tentative: build coarse coordinates") ) {
422 bTransferCoordinates_ =
true;
423 Input(fineLevel,
"Coordinates");
424 }
else if (bTransferCoordinates_) {
425 Input(fineLevel,
"Coordinates");
429 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
430 void TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::Build(Level& fineLevel, Level& coarseLevel)
const {
431 return BuildP(fineLevel, coarseLevel);
434 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
435 void TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::BuildP(Level& fineLevel, Level& coarseLevel)
const {
436 FactoryMonitor m(*
this,
"Build", coarseLevel);
441 const ParameterList& pL = GetParameterList();
442 std::string nspName =
"Nullspace";
443 if(pL.isParameter(
"Nullspace name")) nspName = pL.get<std::string>(
"Nullspace name");
445 auto A = Get< RCP<Matrix> > (fineLevel,
"A");
446 auto aggregates = Get< RCP<Aggregates_kokkos> > (fineLevel,
"Aggregates");
447 auto amalgInfo = Get< RCP<AmalgamationInfo_kokkos> > (fineLevel,
"UnAmalgamationInfo");
448 auto fineNullspace = Get< RCP<MultiVector> > (fineLevel, nspName);
449 auto coarseMap = Get< RCP<const Map> > (fineLevel,
"CoarseMap");
450 RCP<RealValuedMultiVector> fineCoords;
451 if(bTransferCoordinates_) {
452 fineCoords = Get< RCP<RealValuedMultiVector> >(fineLevel,
"Coordinates");
455 RCP<Matrix> Ptentative;
456 RCP<MultiVector> coarseNullspace;
457 RCP<RealValuedMultiVector> coarseCoords;
459 if(bTransferCoordinates_) {
460 ArrayView<const GO> elementAList = coarseMap->getNodeElementList();
461 GO indexBase = coarseMap->getIndexBase();
464 if (rcp_dynamic_cast<const StridedMap>(coarseMap) != Teuchos::null)
465 blkSize = rcp_dynamic_cast<
const StridedMap>(coarseMap)->getFixedBlockSize();
467 Array<GO> elementList;
468 ArrayView<const GO> elementListView;
472 elementListView = elementAList;
475 auto numElements = elementAList.size() / blkSize;
477 elementList.resize(numElements);
480 for (
LO i = 0; i < Teuchos::as<LO>(numElements); i++)
481 elementList[i] = (elementAList[i*blkSize]-indexBase)/blkSize + indexBase;
483 elementListView = elementList;
486 auto uniqueMap = fineCoords->getMap();
488 elementListView, indexBase, coarseMap->getComm());
489 coarseCoords = RealValuedMultiVectorFactory::Build(coarseCoordMap, fineCoords->getNumVectors());
492 RCP<RealValuedMultiVector> ghostedCoords = fineCoords;
493 if (aggregates->AggregatesCrossProcessors()) {
494 auto nonUniqueMap = aggregates->GetMap();
495 auto importer = ImportFactory::Build(uniqueMap, nonUniqueMap);
497 ghostedCoords = RealValuedMultiVectorFactory::Build(nonUniqueMap, fineCoords->getNumVectors());
503 auto aggGraph = aggregates->GetGraph();
504 auto numAggs = aggGraph.numRows();
506 auto fineCoordsView = fineCoords ->template getLocalView<DeviceType>();
507 auto coarseCoordsView = coarseCoords->template getLocalView<DeviceType>();
511 SubFactoryMonitor m2(*
this,
"AverageCoords", coarseLevel);
513 const auto dim = fineCoords->getNumVectors();
515 typename AppendTrait<decltype(fineCoordsView), Kokkos::RandomAccess>::type fineCoordsRandomView = fineCoordsView;
516 for (
size_t j = 0; j < dim; j++) {
518 KOKKOS_LAMBDA(
const LO i) {
522 auto aggregate = aggGraph.rowConst(i);
524 coordinate_type sum = 0.0;
525 for (
size_t colID = 0; colID < static_cast<size_t>(aggregate.length); colID++)
526 sum += fineCoordsRandomView(aggregate(colID),j);
528 coarseCoordsView(i,j) = sum / aggregate.length;
534 if (!aggregates->AggregatesCrossProcessors())
535 BuildPuncoupled(coarseLevel, A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace, coarseLevel.GetLevelID());
537 BuildPcoupled (A, aggregates, amalgInfo, fineNullspace, coarseMap, Ptentative, coarseNullspace);
547 if (A->IsView(
"stridedMaps") ==
true)
548 Ptentative->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), coarseMap);
550 Ptentative->CreateView(
"stridedMaps", Ptentative->getRangeMap(), coarseMap);
552 if(bTransferCoordinates_) {
553 Set(coarseLevel,
"Coordinates", coarseCoords);
555 Set(coarseLevel,
"Nullspace", coarseNullspace);
556 Set(coarseLevel,
"P", Ptentative);
559 RCP<ParameterList> params =
rcp(
new ParameterList());
560 params->set(
"printLoadBalancingInfo",
true);
565 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
566 void TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::
567 BuildPuncoupled(Level& coarseLevel, RCP<Matrix> A, RCP<Aggregates_kokkos> aggregates,
568 RCP<AmalgamationInfo_kokkos> amalgInfo, RCP<MultiVector> fineNullspace,
569 RCP<const Map> coarseMap, RCP<Matrix>& Ptentative,
570 RCP<MultiVector>& coarseNullspace,
const int levelID)
const {
571 auto rowMap = A->getRowMap();
572 auto colMap = A->getColMap();
574 const size_t numRows = rowMap->getNodeNumElements();
575 const size_t NSDim = fineNullspace->getNumVectors();
577 typedef Kokkos::ArithTraits<SC> ATS;
578 using impl_ATS = Kokkos::ArithTraits<typename ATS::val_type>;
579 const SC zero = ATS::zero(), one = ATS::one();
583 typename Aggregates_kokkos::local_graph_type aggGraph;
585 SubFactoryMonitor m2(*
this,
"Get Aggregates graph", coarseLevel);
586 aggGraph = aggregates->GetGraph();
588 auto aggRows = aggGraph.row_map;
589 auto aggCols = aggGraph.entries;
596 SubFactoryMonitor m2(*
this,
"Check good map", coarseLevel);
597 goodMap = isGoodMap(*rowMap, *colMap);
601 "MueLu: TentativePFactory_kokkos: for now works only with good maps "
602 "(i.e. \"matching\" row and column maps)");
611 LO fullBlockSize, blockID, stridingOffset, stridedBlockSize;
613 amalgInfo->GetStridingInformation(fullBlockSize, blockID, stridingOffset, stridedBlockSize, indexBase);
614 GO globalOffset = amalgInfo->GlobalOffset();
617 auto procWinner = aggregates->GetProcWinner() ->template getLocalView<DeviceType>();
618 auto vertex2AggId = aggregates->GetVertex2AggId()->template getLocalView<DeviceType>();
619 const size_t numAggregates = aggregates->GetNumAggregates();
621 int myPID = aggregates->GetMap()->getComm()->getRank();
626 typedef typename Aggregates_kokkos::aggregates_sizes_type::non_const_type AggSizeType;
627 AggSizeType aggDofSizes;
629 if (stridedBlockSize == 1) {
630 SubFactoryMonitor m2(*
this,
"Calc AggSizes", coarseLevel);
633 aggDofSizes = AggSizeType(
"agg_dof_sizes", numAggregates+1);
635 auto sizesConst = aggregates->ComputeAggregateSizes();
639 SubFactoryMonitor m2(*
this,
"Calc AggSizes", coarseLevel);
642 aggDofSizes = AggSizeType(
"agg_dof_sizes", numAggregates + 1);
644 auto nodeMap = aggregates->GetMap()->getLocalMap();
645 auto dofMap = colMap->getLocalMap();
648 KOKKOS_LAMBDA(
const LO agg) {
649 auto aggRowView = aggGraph.rowConst(agg);
652 for (
LO colID = 0; colID < aggRowView.length; colID++) {
653 GO nodeGID = nodeMap.getGlobalElement(aggRowView(colID));
655 for (
LO k = 0; k < stridedBlockSize; k++) {
656 GO dofGID = (nodeGID - indexBase) * fullBlockSize + k + indexBase + globalOffset + stridingOffset;
658 if (dofMap.getLocalElement(dofGID) != INVALID)
662 aggDofSizes(agg+1) = size;
669 ReduceMaxFunctor<LO,decltype(aggDofSizes)> reduceMax(aggDofSizes);
670 Kokkos::parallel_reduce(
"MueLu:TentativePF:Build:max_agg_size", range_type(0, aggDofSizes.extent(0)), reduceMax, maxAggSize);
674 Kokkos::parallel_scan(
"MueLu:TentativePF:Build:aggregate_sizes:stage1_scan", range_type(0,numAggregates+1),
675 KOKKOS_LAMBDA(
const LO i,
LO& update,
const bool& final_pass) {
676 update += aggDofSizes(i);
678 aggDofSizes(i) = update;
685 SubFactoryMonitor m2(*
this,
"Create Agg2RowMap", coarseLevel);
687 AggSizeType aggOffsets(Kokkos::ViewAllocateWithoutInitializing(
"aggOffsets"), numAggregates);
690 Kokkos::parallel_for(
"MueLu:TentativePF:Build:createAgg2RowMap", range_type(0, vertex2AggId.extent(0)),
691 KOKKOS_LAMBDA(
const LO lnode) {
692 if (procWinner(lnode, 0) == myPID) {
694 auto aggID = vertex2AggId(lnode,0);
696 auto offset = Kokkos::atomic_fetch_add( &aggOffsets(aggID), stridedBlockSize );
700 for (
LO k = 0; k < stridedBlockSize; k++)
701 agg2RowMapLO(offset + k) = lnode*stridedBlockSize + k;
708 coarseNullspace = MultiVectorFactory::Build(coarseMap, NSDim);
711 auto fineNS = fineNullspace ->template getLocalView<DeviceType>();
712 auto coarseNS = coarseNullspace->template getLocalView<DeviceType>();
717 typedef typename local_matrix_type::row_map_type::non_const_type rows_type;
718 typedef typename local_matrix_type::index_type::non_const_type cols_type;
719 typedef typename local_matrix_type::values_type::non_const_type vals_type;
724 status_type status(
"status");
726 typename AppendTrait<decltype(fineNS), Kokkos::RandomAccess>::type fineNSRandom = fineNS;
727 typename AppendTrait<status_type, Kokkos::Atomic> ::type statusAtomic = status;
733 const ParameterList& pL = GetParameterList();
734 const bool& doQRStep = pL.get<
bool>(
"tentative: calculate qr");
736 GetOStream(
Runtime1) <<
"TentativePFactory : bypassing local QR phase" << std::endl;
738 GetOStream(
Warnings0) <<
"TentativePFactor : for nontrivial nullspace, this may degrade performance" << std::endl;
747 SubFactoryMonitor m2(*
this,
"Stage 1 (LocalQR)", coarseLevel);
752 rows = rows_type(
"Ptent_rows", numRows+1);
753 cols = cols_type(Kokkos::ViewAllocateWithoutInitializing(
"Ptent_cols"), numRows);
754 vals = vals_type(Kokkos::ViewAllocateWithoutInitializing(
"Ptent_vals"), numRows);
765 auto agg = thread.league_rank();
768 LO aggSize = aggRows(agg+1) - aggRows(agg);
773 auto norm = impl_ATS::magnitude(zero);
778 for (decltype(aggSize) k = 0; k < aggSize; k++) {
779 auto dnorm = impl_ATS::magnitude(fineNSRandom(agg2RowMapLO(aggRows(agg)+k),0));
786 statusAtomic(1) =
true;
791 coarseNS(agg, 0) = norm;
794 for (decltype(aggSize) k = 0; k < aggSize; k++) {
795 LO localRow = agg2RowMapLO(aggRows(agg)+k);
796 SC localVal = fineNSRandom(agg2RowMapLO(aggRows(agg)+k),0) / norm;
798 rows(localRow+1) = localRow+1;
799 cols(localRow) = agg;
800 vals(localRow) = localVal;
805 typename status_type::HostMirror statusHost = Kokkos::create_mirror_view(status);
807 for (decltype(statusHost.size()) i = 0; i < statusHost.size(); i++)
809 std::ostringstream oss;
810 oss <<
"MueLu::TentativePFactory::MakeTentative: ";
812 case 0: oss <<
"!goodMap is not implemented";
break;
813 case 1: oss <<
"fine level NS part has a zero column";
break;
815 throw Exceptions::RuntimeError(oss.str());
821 auto agg = thread.league_rank();
824 LO aggSize = aggRows(agg+1) - aggRows(agg);
827 coarseNS(agg, 0) = one;
830 for (decltype(aggSize) k = 0; k < aggSize; k++) {
831 LO localRow = agg2RowMapLO(aggRows(agg)+k);
832 SC localVal = fineNSRandom(agg2RowMapLO(aggRows(agg)+k),0);
834 rows(localRow+1) = localRow+1;
835 cols(localRow) = agg;
836 vals(localRow) = localVal;
852 size_t nnzEstimate = numRows * NSDim;
853 rows_type rowsAux(
"Ptent_aux_rows", numRows+1);
854 cols_type colsAux(
"Ptent_aux_cols", nnzEstimate);
855 vals_type valsAux(
"Ptent_aux_vals", nnzEstimate);
856 rows = rows_type(
"Ptent_rows", numRows+1);
859 SubFactoryMonitor m2(*
this,
"Stage 0 (InitViews)", coarseLevel);
864 KOKKOS_LAMBDA(
const LO row) {
865 rowsAux(row) = row*NSDim;
868 KOKKOS_LAMBDA(
const LO j) {
869 colsAux(j) = INVALID;
870 valsAux(j) = impl_ATS::zero();
875 SubFactoryMonitor m2 = SubFactoryMonitor(*
this, doQRStep ?
"Stage 1 (LocalQR)" :
"Stage 1 (Fill coarse nullspace and tentative P)", coarseLevel);
881 decltype(aggDofSizes ), decltype(maxAggSize), decltype(agg2RowMapLO),
882 decltype(statusAtomic), decltype(rows), decltype(rowsAux), decltype(colsAux),
884 localQRFunctor(fineNSRandom, coarseNS, aggDofSizes, maxAggSize, agg2RowMapLO, statusAtomic,
885 rows, rowsAux, colsAux, valsAux, doQRStep);
889 typename status_type::HostMirror statusHost = Kokkos::create_mirror_view(status);
891 for (decltype(statusHost.size()) i = 0; i < statusHost.size(); i++)
893 std::ostringstream oss;
894 oss <<
"MueLu::TentativePFactory::MakeTentative: ";
896 case 0: oss <<
"!goodMap is not implemented";
break;
897 case 1: oss <<
"fine level NS part has a zero column";
break;
899 throw Exceptions::RuntimeError(oss.str());
906 cols = decltype(cols)(
"Ptent_cols", nnz);
907 vals = decltype(vals)(
"Ptent_vals", nnz);
910 SubFactoryMonitor m2(*
this,
"Stage 2 (CompressRows)", coarseLevel);
912 Kokkos::parallel_scan(
"MueLu:TentativePF:Build:compress_rows", range_type(0,numRows+1),
913 KOKKOS_LAMBDA(
const LO i,
LO& upd,
const bool&
final) {
921 SubFactoryMonitor m2(*
this,
"Stage 2 (CompressCols)", coarseLevel);
927 KOKKOS_LAMBDA(
const LO i) {
928 LO rowStart = rows(i);
931 for (
auto j = rowsAux(i); j < rowsAux(i+1); j++)
932 if (colsAux(j) != INVALID) {
933 cols(rowStart+lnnz) = colsAux(j);
934 vals(rowStart+lnnz) = valsAux(j);
941 GetOStream(
Runtime1) <<
"TentativePFactory : aggregates do not cross process boundaries" << std::endl;
945 SubFactoryMonitor m2(*
this,
"Stage 3 (LocalMatrix+FillComplete)", coarseLevel);
947 local_matrix_type lclMatrix = local_matrix_type(
"A", numRows, coarseMap->getNodeNumElements(), nnz, vals, rows, cols);
950 RCP<ParameterList> FCparams;
951 if (pL.isSublist(
"matrixmatrix: kernel params"))
952 FCparams =
rcp(
new ParameterList(pL.sublist(
"matrixmatrix: kernel params")));
954 FCparams =
rcp(
new ParameterList);
957 FCparams->set(
"compute global constants", FCparams->get(
"compute global constants",
false));
958 FCparams->set(
"Timer Label", std::string(
"MueLu::TentativeP-") +
toString(levelID));
960 auto PtentCrs = CrsMatrixFactory::Build(lclMatrix, rowMap, coarseMap, coarseMap, A->getDomainMap());
961 Ptentative =
rcp(
new CrsMatrixWrap(PtentCrs));
965 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
966 void TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::
967 BuildPcoupled(RCP<Matrix> , RCP<Aggregates_kokkos> ,
968 RCP<AmalgamationInfo_kokkos> , RCP<MultiVector> ,
969 RCP<const Map> , RCP<Matrix>& ,
970 RCP<MultiVector>& )
const {
971 throw Exceptions::RuntimeError(
"MueLu: Construction of coupled tentative P is not implemented");
974 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class DeviceType>
975 bool TentativePFactory_kokkos<Scalar,LocalOrdinal,GlobalOrdinal,Kokkos::Compat::KokkosDeviceWrapperNode<DeviceType>>::
976 isGoodMap(
const Map& rowMap,
const Map& colMap)
const {
977 auto rowLocalMap = rowMap.getLocalMap();
978 auto colLocalMap = colMap.getLocalMap();
980 const size_t numRows = rowLocalMap.getNodeNumElements();
981 const size_t numCols = colLocalMap.getNodeNumElements();
983 if (numCols < numRows)
988 KOKKOS_LAMBDA(
const LO i,
size_t &diff) {
989 diff += (rowLocalMap.getGlobalElement(i) != colLocalMap.getGlobalElement(i));
992 return (numDiff == 0);
997 #define MUELU_TENTATIVEPFACTORY_KOKKOS_SHORT
998 #endif // HAVE_MUELU_KOKKOS_REFACTOR
999 #endif // MUELU_TENTATIVEPFACTORY_KOKKOS_DEF_HPP
Important warning messages (one line)
MueLu::DefaultLocalOrdinal LocalOrdinal
void parallel_for(const ExecPolicy &policy, const FunctorType &functor, const std::string &str="", typename Impl::enable_if< Kokkos::Impl::is_execution_policy< ExecPolicy >::value >::type *=0)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
void parallel_reduce(const std::string &label, const PolicyType &policy, const FunctorType &functor, ReturnType &return_value, typename Impl::enable_if< Kokkos::Impl::is_execution_policy< PolicyType >::value >::type *=0)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
KOKKOS_INLINE_FUNCTION Kokkos::complex< RealType > pow(const complex< RealType > &x, const RealType &e)
static const NoFactory * get()
void deep_copy(const View< DT, DP...> &dst, typename ViewTraits< DT, DP...>::const_value_type &value, typename std::enable_if< std::is_same< typename ViewTraits< DT, DP...>::specialize, void >::value >::type *=0)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
KOKKOS_INLINE_FUNCTION Kokkos::complex< RealType > sqrt(const complex< RealType > &x)
KOKKOS_FORCEINLINE_FUNCTION constexpr pair< T1, T2 > make_pair(T1 x, T2 y)
Description of what is happening (more verbose)
#define SET_VALID_ENTRY(name)