46 #ifndef MUELU_AMALGAMATIONFACTORY_KOKKOS_DEF_HPP
47 #define MUELU_AMALGAMATIONFACTORY_KOKKOS_DEF_HPP
49 #include <Xpetra_Matrix.hpp>
55 #ifdef HAVE_MUELU_KOKKOS_REFACTOR
57 #include "MueLu_AmalgamationInfo_kokkos.hpp"
61 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
62 RCP<const ParameterList> AmalgamationFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
63 GetValidParameterList()
const {
64 RCP<ParameterList> validParamList =
rcp(
new ParameterList());
65 validParamList->set< RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory of the matrix A");
66 return validParamList;
69 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 void AmalgamationFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
71 DeclareInput(Level ¤tLevel)
const {
72 Input(currentLevel,
"A");
75 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 void AmalgamationFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
77 Build(Level ¤tLevel)
const
79 FactoryMonitor m(*
this,
"Build", currentLevel);
81 RCP<Matrix> A = Get< RCP<Matrix> >(currentLevel,
"A");
86 LO nStridedOffset = 0;
87 LO stridedblocksize = fullblocksize;
92 if (A->IsView(
"stridedMaps") && Teuchos::rcp_dynamic_cast<
const StridedMap>(A->getRowMap(
"stridedMaps")) != Teuchos::null) {
93 Xpetra::viewLabel_t oldView = A->SwitchToView(
"stridedMaps");
94 RCP<const StridedMap> stridedRowMap = Teuchos::rcp_dynamic_cast<
const StridedMap>(A->getRowMap());
95 TEUCHOS_TEST_FOR_EXCEPTION(stridedRowMap == Teuchos::null,Exceptions::BadCast,
"MueLu::CoalesceFactory::Build: cast to strided row map failed.");
96 fullblocksize = stridedRowMap->getFixedBlockSize();
97 offset = stridedRowMap->getOffset();
98 blockid = stridedRowMap->getStridedBlockId();
101 std::vector<size_t> stridingInfo = stridedRowMap->getStridingData();
102 for (
size_t j = 0; j < Teuchos::as<size_t>(blockid); j++)
103 nStridedOffset += stridingInfo[j];
104 stridedblocksize = Teuchos::as<LocalOrdinal>(stridingInfo[blockid]);
107 stridedblocksize = fullblocksize;
109 oldView = A->SwitchToView(oldView);
110 GetOStream(
Runtime1) <<
"AmalagamationFactory::Build():" <<
" found fullblocksize=" << fullblocksize <<
" and stridedblocksize=" << stridedblocksize <<
" from strided maps. offset=" << offset << std::endl;
113 GetOStream(
Warnings0) <<
"AmalagamationFactory::Build(): no striding information available. Use blockdim=1 with offset=0" << std::endl;
120 RCP<const Map> uniqueMap, nonUniqueMap;
121 RCP<AmalgamationInfo_kokkos> amalgamationData;
122 RCP<Array<LO> > rowTranslation = Teuchos::null;
123 RCP<Array<LO> > colTranslation = Teuchos::null;
125 if (fullblocksize > 1) {
131 RCP<Array<LO> > theRowTranslation =
rcp(
new Array<LO>);
132 RCP<Array<LO> > theColTranslation =
rcp(
new Array<LO>);
133 AmalgamateMap(*(A->getRowMap()), *A, uniqueMap, *theRowTranslation);
134 AmalgamateMap(*(A->getColMap()), *A, nonUniqueMap, *theColTranslation);
136 amalgamationData =
rcp(
new AmalgamationInfo_kokkos(theRowTranslation,
147 amalgamationData =
rcp(
new AmalgamationInfo_kokkos(rowTranslation,
160 Set(currentLevel,
"UnAmalgamationInfo", amalgamationData);
163 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
164 void AmalgamationFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
165 AmalgamateMap(
const Map& sourceMap,
const Matrix& A, RCP<const Map>& amalgamatedMap, Array<LO>& translation) {
166 typedef typename ArrayView<const GO>::size_type size_type;
167 typedef std::map<GO,size_type> container;
169 GO indexBase = sourceMap.getIndexBase();
170 ArrayView<const GO> elementAList = sourceMap.getNodeElementList();
171 size_type numElements = elementAList.size();
175 LO blkSize = A.GetFixedBlockSize();
176 if (A.IsView(
"stridedMaps") ==
true) {
180 offset = strMap->getOffset();
181 blkSize = Teuchos::as<const LO>(strMap->getFixedBlockSize());
184 Array<GO> elementList(numElements);
185 translation.resize(numElements);
187 size_type numRows = 0;
188 for (size_type
id = 0;
id < numElements;
id++) {
189 GO dofID = elementAList[id];
190 GO nodeID = AmalgamationFactory_kokkos::DOFGid2NodeId(dofID, blkSize, offset, indexBase);
192 typename container::iterator it = filter.find(nodeID);
193 if (it == filter.end()) {
194 filter[nodeID] = numRows;
196 translation[id] = numRows;
197 elementList[numRows] = nodeID;
202 translation[id] = it->second;
205 elementList.resize(numRows);
211 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
212 const GlobalOrdinal AmalgamationFactory_kokkos<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
216 return globalblockid;
221 #endif // HAVE_MUELU_KOKKOS_REFACTOR
Important warning messages (one line)
MueLu::DefaultLocalOrdinal LocalOrdinal
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Description of what is happening (more verbose)