10 #ifndef MUELU_AMALGAMATIONFACTORY_DEF_HPP
11 #define MUELU_AMALGAMATIONFACTORY_DEF_HPP
18 #include "MueLu_AmalgamationInfo.hpp"
23 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
26 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
29 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
33 return validParamList;
36 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
38 Input(currentLevel,
"A");
41 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
45 RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel,
"A");
67 LO nStridedOffset = 0;
68 LO stridedblocksize = fullblocksize;
69 LO storageblocksize = A->GetStorageBlockSize();
74 if (A->IsView(
"stridedMaps") && Teuchos::rcp_dynamic_cast<
const StridedMap>(A->getRowMap(
"stridedMaps")) != Teuchos::null) {
78 fullblocksize = stridedRowMap->getFixedBlockSize();
79 offset = DOFGidOffset(stridedRowMap);
80 blockid = stridedRowMap->getStridedBlockId();
83 std::vector<size_t> stridingInfo = stridedRowMap->getStridingData();
84 for (
size_t j = 0; j < Teuchos::as<size_t>(blockid); j++)
85 nStridedOffset += stridingInfo[j];
86 stridedblocksize = Teuchos::as<LocalOrdinal>(stridingInfo[blockid]);
89 stridedblocksize = fullblocksize;
93 TEUCHOS_TEST_FOR_EXCEPTION(fullblocksize % storageblocksize != 0,
Exceptions::RuntimeError,
"AmalgamationFactory: fullblocksize needs to be a multiple of A->GetStorageBlockSize()");
94 fullblocksize /= storageblocksize;
95 stridedblocksize /= storageblocksize;
97 oldView = A->SwitchToView(oldView);
98 GetOStream(
Runtime1) <<
"AmalagamationFactory::Build():"
99 <<
" found fullblocksize=" << fullblocksize <<
" and stridedblocksize=" << stridedblocksize <<
" from strided maps. offset=" << offset << std::endl;
102 GetOStream(
Warnings0) <<
"AmalagamationFactory::Build(): no striding information available. Use blockdim=1 with offset=0" << std::endl;
114 if (fullblocksize > 1) {
122 AmalgamateMap(*(A->getRowMap()), *A, uniqueMap, *theRowTranslation);
123 AmalgamateMap(*(A->getColMap()), *A, nonUniqueMap, *theColTranslation);
149 Set(currentLevel,
"UnAmalgamationInfo", amalgamationData);
152 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
155 typedef std::unordered_map<GO, size_type> container;
157 GO indexBase = sourceMap.getIndexBase();
159 size_type numElements = elementAList.
size();
163 LO blkSize = A.GetFixedBlockSize() / A.GetStorageBlockSize();
164 if (A.IsView(
"stridedMaps") ==
true) {
168 offset = DOFGidOffset(strMap);
169 blkSize = Teuchos::as<const LO>(strMap->getFixedBlockSize());
173 translation.
resize(numElements);
175 size_type numRows = 0;
176 for (size_type
id = 0;
id < numElements;
id++) {
177 GO dofID = elementAList[id];
180 typename container::iterator it = filter.find(nodeID);
181 if (it == filter.end()) {
182 filter[nodeID] = numRows;
184 translation[id] = numRows;
185 elementList[numRows] = nodeID;
190 translation[id] = it->second;
193 elementList.
resize(numRows);
198 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
202 return globalblockid;
205 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
207 GlobalOrdinal offset = stridedRowMap->getMinAllGlobalIndex();
209 size_t nStridedOffset = 0;
210 for (
int j = 0; j < stridedRowMap->getStridedBlockId(); j++) {
211 nStridedOffset += stridedRowMap->getStridingData()[j];
213 const GO goStridedOffset = Teuchos::as<const GO>(nStridedOffset);
215 offset -= goStridedOffset + stridedRowMap->getIndexBase();
Important warning messages (one line)
void DeclareInput(Level ¤tLevel) const override
Input.
Exception indicating invalid cast attempted.
MueLu::DefaultLocalOrdinal LocalOrdinal
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static const GlobalOrdinal DOFGidOffset(RCP< const StridedMap > stridedMap)
Method to calculate the global (row) id offset from scratch.
static const GlobalOrdinal DOFGid2NodeId(GlobalOrdinal gid, LocalOrdinal blockSize, const GlobalOrdinal offset, const GlobalOrdinal indexBase)
Translate global (row/column) id to global amalgamation block id.
AmalgamationFactory()
Constructor.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
static void AmalgamateMap(const Map &sourceMap, const Matrix &A, RCP< const Map > &amalgamatedMap, Array< LO > &translation)
Method to create merged map for systems of PDEs.
void resize(size_type new_size, const value_type &x=value_type())
void Build(Level ¤tLevel) const override
Build an object with this factory.
virtual ~AmalgamationFactory()
Destructor.
RCP< const ParameterList > GetValidParameterList() const override
Return a const parameter list of valid parameters that setParameterList() will accept.
Exception throws to report errors in the internal logical of the program.
Description of what is happening (more verbose)
minimal container class for storing amalgamation information