MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_ReorderBlockAFactory_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_REORDERBLOCKAFACTORY_DEF_HPP_
11#define MUELU_REORDERBLOCKAFACTORY_DEF_HPP_
12
14
15#include <Xpetra_BlockReorderManager.hpp>
16#include <Xpetra_BlockedCrsMatrix.hpp>
17#include <Xpetra_ReorderedBlockedCrsMatrix.hpp>
18#include <Xpetra_MapExtractor.hpp>
19#include <Xpetra_MapExtractorFactory.hpp>
20#include <Xpetra_Matrix.hpp>
21#include <Xpetra_MatrixUtils.hpp>
22
23#include "MueLu_Level.hpp"
24#include "MueLu_Monitor.hpp"
25
26namespace MueLu {
27
28template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30 RCP<ParameterList> validParamList = rcp(new ParameterList());
31
32 validParamList->set<RCP<const FactoryBase> >("A", MueLu::NoFactory::getRCP(), "Generating factory for A.");
33
34 validParamList->set<std::string>("Reorder Type", "", "String describing the reordering of blocks");
35
36 // TODO not very elegant.
37 validParamList->set<RCP<const FactoryBase> >("Map1", Teuchos::null, "Generating factory of the fine level map associated with the (1,1) block in your n x n block matrix.");
38 validParamList->set<RCP<const FactoryBase> >("Map2", Teuchos::null, "Generating factory of the fine level map associated with the (2,2) block in your n x n block matrix.");
39 validParamList->set<RCP<const FactoryBase> >("Map3", Teuchos::null, "Generating factory of the fine level map associated with the (3,3) block in your n x n block matrix.");
40 validParamList->set<RCP<const FactoryBase> >("Map4", Teuchos::null, "Generating factory of the fine level map associated with the (4,4) block in your n x n block matrix.");
41 validParamList->set<RCP<const FactoryBase> >("Map5", Teuchos::null, "Generating factory of the fine level map associated with the (5,5) block in your n x n block matrix.");
42 validParamList->set<RCP<const FactoryBase> >("Map6", Teuchos::null, "Generating factory of the fine level map associated with the (6,6) block in your n x n block matrix.");
43 validParamList->set<RCP<const FactoryBase> >("Map7", Teuchos::null, "Generating factory of the fine level map associated with the (7,7) block in your n x n block matrix.");
44 validParamList->set<RCP<const FactoryBase> >("Map8", Teuchos::null, "Generating factory of the fine level map associated with the (8,8) block in your n x n block matrix.");
45 validParamList->set<RCP<const FactoryBase> >("Map9", Teuchos::null, "Generating factory of the fine level map associated with the (9,9) block in your n x n block matrix.");
46
47 return validParamList;
48}
49
50template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
52 Input(currentLevel, "A");
53}
54
55template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57 FactoryMonitor m(*this, "ReorderBlockA factory", currentLevel);
58
59 const ParameterList& pL = GetParameterList();
60 std::string reorderStr = pL.get<std::string>("Reorder Type");
61
62 RCP<Matrix> Ain = Get<RCP<Matrix> >(currentLevel, "A");
63
64 RCP<BlockedCrsMatrix> A = rcp_dynamic_cast<BlockedCrsMatrix>(Ain);
65
66 // special case: we get a single block CrsMatrix object on the finest level and
67 // split it into a nxn blocked operator
68 if (A == Teuchos::null && currentLevel.GetLevelID() == 0) {
69 GetOStream(Warnings0) << "Split input matrix (Warning: this is a rather expensive operation)" << std::endl;
70
71 std::vector<Teuchos::RCP<const Map> > xmaps;
72
73 for (int it = 1; it < 10; it++) {
74 std::stringstream ss;
75 ss << "Map" << it;
76 if (currentLevel.IsAvailable(ss.str(), NoFactory::get())) {
77 RCP<const Map> submap = currentLevel.Get<RCP<const Map> >(ss.str(), NoFactory::get());
78 GetOStream(Runtime1) << "Use user-given submap #" << it << ": length dimension=" << submap->getGlobalNumElements() << std::endl;
79 xmaps.push_back(submap);
80 }
81 }
82
83 bool bThyraMode = false; // no support for Thyra mode (yet)
84 RCP<const MapExtractor> map_extractor = Xpetra::MapExtractorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(Ain->getRowMap(), xmaps, bThyraMode);
85
86 // split null space vectors
87 // TODO: if he matrix blocks have different striding, this could be quite complicated
88 // RCP<MultiVector> nullspace1 = map_extractor->ExtractVector(nullspace,0);
89 // RCP<MultiVector> nullspace2 = map_extractor->ExtractVector(nullspace,1);
90
91 Teuchos::RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > bOp =
92 Xpetra::MatrixUtils<Scalar, LocalOrdinal, GlobalOrdinal, Node>::SplitMatrix(*Ain, map_extractor, map_extractor, Teuchos::null, bThyraMode);
93
94 TEUCHOS_TEST_FOR_EXCEPTION(Ain->getGlobalNumRows() != bOp->getGlobalNumRows(), Exceptions::RuntimeError, "Split operator not consistent with input operator (different number of rows).");
95 TEUCHOS_TEST_FOR_EXCEPTION(Ain->getLocalNumRows() != bOp->getLocalNumRows(), Exceptions::RuntimeError, "Split operator not consistent with input operator (different number of node rows).");
96 TEUCHOS_TEST_FOR_EXCEPTION(Ain->getLocalNumEntries() != bOp->getLocalNumEntries(), Exceptions::RuntimeError, "Split operator not consistent with input operator (different number of local entries).");
97 TEUCHOS_TEST_FOR_EXCEPTION(Ain->getGlobalNumEntries() != bOp->getGlobalNumEntries(), Exceptions::RuntimeError, "Split operator not consistent with input operator (different number of global entries).");
98
99 A = bOp;
100 }
101
102 // we have a blocked operator as input
103 TEUCHOS_TEST_FOR_EXCEPTION(A.is_null(), Exceptions::BadCast, "Input matrix A is not a BlockedCrsMatrix.");
104 GetOStream(Statistics1) << "Got a " << A->Rows() << "x" << A->Cols() << " blocked operator as input" << std::endl;
105
106 // if we have a blocked operator and a reordering string, create a nested blocked operator, if not skip the process
107 if (reorderStr.empty()) {
108 GetOStream(Statistics1) << "No reordering information provided. Skipping reordering of A." << std::endl;
109 } else {
110 Teuchos::RCP<const Xpetra::BlockReorderManager> brm = Xpetra::blockedReorderFromString(reorderStr);
111 GetOStream(Debug) << "Reordering A using " << brm->toString() << std::endl;
112
113 Teuchos::RCP<const ReorderedBlockedCrsMatrix> brop =
114 Teuchos::rcp_dynamic_cast<const ReorderedBlockedCrsMatrix>(
115 Xpetra::buildReorderedBlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(brm, A));
116
117 TEUCHOS_TEST_FOR_EXCEPTION(brop.is_null(), Exceptions::RuntimeError,
118 "Block reordering of " << A->Rows() << "x" << A->Cols()
119 << " blocked operator failed.");
120
121 GetOStream(Statistics1) << "Reordering A using " << brm->toString() << " block gives a " << brop->Rows() << "x"
122 << brop->Cols() << " blocked operator" << std::endl;
123 GetOStream(Debug) << "Reordered operator has " << brop->getRangeMap()->getGlobalNumElements() << " rows and "
124 << brop->getDomainMap()->getGlobalNumElements() << " columns" << std::endl;
125 GetOStream(Debug) << "Reordered operator: Use of Thyra style gids = "
126 << brop->getRangeMapExtractor()->getThyraMode() << std::endl;
127
128 // get rid of const (we expect non-const operators stored in Level)
129 Teuchos::RCP<ReorderedBlockedCrsMatrix> bret =
130 Teuchos::rcp_const_cast<ReorderedBlockedCrsMatrix>(brop);
131
132 A = bret;
133 }
134
135 currentLevel.Set("A", Teuchos::rcp_dynamic_cast<Matrix>(A), this);
136}
137
138} // namespace MueLu
139
140#endif /* MUELU_REORDERBLOCKAFACTORY_DEF_HPP_ */
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
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)....
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
static const RCP< const NoFactory > getRCP()
Static Get() functions.
static const NoFactory * get()
RCP< const ParameterList > GetValidParameterList() const
Input.
void Build(Level &currentLevel) const
Build an object with this factory.
void DeclareInput(Level &currentLevel) const
Specifies the data that this class needs, and the factories that generate that data.
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Debug
Print additional debugging information.
@ Statistics1
Print more statistics.
@ Runtime1
Description of what is happening (more verbose)