MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_EdgeProlongatorPatternFactory_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_EDGEPROLONGATORPATTERNFACTORY_DEF_HPP
11#define MUELU_EDGEPROLONGATORPATTERNFACTORY_DEF_HPP
12
13#include <Xpetra_Matrix.hpp>
14#include <Xpetra_MatrixMatrix.hpp>
15
17
18#include "MueLu_Monitor.hpp"
19#include "Teuchos_ScalarTraits.hpp"
20
21namespace MueLu {
22
23template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
25 RCP<ParameterList> validParamList = rcp(new ParameterList());
26
27 validParamList->set<RCP<const FactoryBase> >("FineD0", Teuchos::null, "Generating factory for the fine discrete gradient");
28 validParamList->set<RCP<const FactoryBase> >("CoarseD0", Teuchos::null, "Generating factory for the coarse discrete gradient");
29 validParamList->set<RCP<const FactoryBase> >("PnodalEmin", Teuchos::null, "Generating factory for the nodal prolongator");
30
31 return validParamList;
32}
33
34template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
36 Input(fineLevel, "D0", "FineD0");
37 Input(coarseLevel, "D0", "CoarseD0");
38 Input(coarseLevel, "PnodalEmin");
39}
40
41template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43 FactoryMonitor m(*this, "EdgeProlongatorPattern", coarseLevel);
44
45 auto D = Get<RCP<Matrix> >(fineLevel, "D0", "FineD0");
46 auto Dc = Get<RCP<Matrix> >(coarseLevel, "D0", "CoarseD0");
47 auto Pn = Get<RCP<Matrix> >(coarseLevel, "PnodalEmin");
48
49 const auto one = Teuchos::ScalarTraits<Scalar>::one();
50
51 // |FineD| * |Pnodal| * |CoarseD^T|
52
53 RCP<Matrix> absD_absPn_absDcT;
54 RCP<Matrix> absDc;
55 {
56 SubFactoryMonitor m2(*this, "Matrix manipulations", coarseLevel);
57
58 auto absD = MatrixFactory::BuildCopy(D);
59 absD->setAllToScalar(one);
60
61 auto absPn = MatrixFactory::BuildCopy(Pn);
62 absPn->setAllToScalar(one);
63
64 RCP<Matrix> absD_absPn = MatrixMatrix::Multiply(*absD, false, *absPn, false, GetOStream(Statistics2), true, true);
65 absD_absPn->setAllToScalar(one);
66
67 // If we rebalanced then Dc lives on a smaller communicator than D.
68 // Since we need to perform matrix-matrix multiplications with Dc, we construct a version of it that lives on the same communicator.
69 auto comm = absD_absPn->getRowMap()->getComm();
70 if (Dc.is_null() || Dc->getRowMap()->getComm()->getSize() < comm->getSize()) {
71 auto lib = absD_absPn->getRowMap()->lib();
72 if (Dc.is_null()) {
73 Kokkos::View<GlobalOrdinal*, typename Node::memory_space> dummy("", 0);
74 auto big_coarse_nodal_map = MapFactory::Build(lib, -1, dummy, 0, comm);
75 auto big_coarse_edge_map = MapFactory::Build(lib, -1, dummy, 0, comm);
76 auto big_coarse_nodal_colmap = MapFactory::Build(lib, -1, dummy, 0, comm);
77
78 typename Matrix::local_matrix_device_type dummyLocalMatrix;
79 Dc = MatrixFactory::Build(dummyLocalMatrix, big_coarse_edge_map, big_coarse_nodal_colmap, big_coarse_nodal_map, big_coarse_edge_map);
80
81 } else {
82 auto big_coarse_nodal_map = MapFactory::Build(lib, -1, Dc->getDomainMap()->getMyGlobalIndicesDevice(), 0, comm);
83 auto big_coarse_edge_map = MapFactory::Build(lib, -1, Dc->getRangeMap()->getMyGlobalIndicesDevice(), 0, comm);
84 auto big_coarse_nodal_colmap = MapFactory::Build(lib, -1, Dc->getColMap()->getMyGlobalIndicesDevice(), 0, comm);
85
86 Dc = MatrixFactory::Build(Dc->getLocalMatrixDevice(), big_coarse_edge_map, big_coarse_nodal_colmap, big_coarse_nodal_map, big_coarse_edge_map);
87 }
88 }
89 absDc = MatrixFactory::BuildCopy(Dc);
90 absDc->setAllToScalar(one);
91
92 absD_absPn_absDcT = MatrixMatrix::Multiply(*absD_absPn, false, *absDc, true, GetOStream(Statistics2), true, true);
93 }
94
95 RCP<Matrix> filtered;
96 {
97 SubFactoryMonitor m2(*this, "Filtering", coarseLevel);
98 using ATS = KokkosKernels::ArithTraits<typename Matrix::impl_scalar_type>;
99 using magnitudeType = typename ATS::magnitudeType;
100 using magATS = KokkosKernels::ArithTraits<magnitudeType>;
101 auto eps = magATS::epsilon();
102
103 RCP<MultiVector> oneVec = MultiVectorFactory::Build(absDc->getDomainMap(), 1);
104 oneVec->putScalar(one);
105 RCP<MultiVector> singleParent = MultiVectorFactory::Build(absDc->getRowMap(), 1);
106 absDc->apply(*oneVec, *singleParent, Teuchos::NO_TRANS);
107 // ghost singleParent
108 RCP<MultiVector> singleParentGhosted;
109 auto importer = absD_absPn_absDcT->getCrsGraph()->getImporter();
110 if (importer.is_null()) {
111 singleParentGhosted = singleParent;
112 } else {
113 singleParentGhosted = MultiVectorFactory::Build(importer->getTargetMap(), 1);
114 singleParentGhosted->doImport(*singleParent, *importer, Xpetra::INSERT);
115 }
116
117 auto lclSingleParent = singleParentGhosted->getLocalViewDevice(Tpetra::Access::ReadOnly);
118
119 // Filter matrix using criterion
120 filtered = Xpetra::applyFilter_LID(
121 absD_absPn_absDcT,
122 KOKKOS_LAMBDA(const LocalOrdinal row,
123 const LocalOrdinal col,
124 const typename Matrix::impl_scalar_type val) {
125 return ((ATS::magnitude(val - 2.0) < eps) || ((lclSingleParent(col, 0) == 1.0) && (ATS::magnitude(val - 1.0) < eps)));
126 });
127
128 if (IsPrint(Statistics1)) {
129 auto numEntriesBeforeFiltering = absD_absPn_absDcT->getGlobalNumEntries();
130 auto numEntriesAfterFiltering = filtered->getGlobalNumEntries();
131 GetOStream(Statistics1) << "Number of kept entries in filtered pattern for P: " << numEntriesAfterFiltering << "/" << numEntriesBeforeFiltering << std::endl;
132 }
133 }
134
135 auto Ppattern = filtered->getCrsGraph();
136
137 Set(coarseLevel, "Ppattern", Ppattern);
138}
139
140} // namespace MueLu
141
142#endif // MUELU_EDGEPROLONGATORPATTERNFACTORY_DEF_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Statistics1
Print more statistics.