10 #ifndef MUELU_REORDERBLOCKAFACTORY_DEF_HPP_
11 #define MUELU_REORDERBLOCKAFACTORY_DEF_HPP_
18 #include <Xpetra_MapExtractor.hpp>
28 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
34 validParamList->
set<std::string>(
"Reorder Type",
"",
"String describing the reordering of blocks");
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.");
47 return validParamList;
50 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
52 Input(currentLevel,
"A");
55 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
60 std::string reorderStr = pL.
get<std::string>(
"Reorder Type");
62 RCP<Matrix> Ain = Get<RCP<Matrix> >(currentLevel,
"A");
68 if (A == Teuchos::null && currentLevel.GetLevelID() == 0) {
69 GetOStream(
Warnings0) <<
"Split input matrix (Warning: this is a rather expensive operation)" << std::endl;
71 std::vector<Teuchos::RCP<const Map> > xmaps;
73 for (
int it = 1; it < 10; it++) {
78 GetOStream(
Runtime1) <<
"Use user-given submap #" << it <<
": length dimension=" << submap->getGlobalNumElements() << std::endl;
79 xmaps.push_back(submap);
83 bool bThyraMode =
false;
104 GetOStream(
Statistics1) <<
"Got a " << A->Rows() <<
"x" << A->Cols() <<
" blocked operator as input" << std::endl;
107 if (reorderStr.empty()) {
108 GetOStream(
Statistics1) <<
"No reordering information provided. Skipping reordering of A." << std::endl;
111 GetOStream(
Debug) <<
"Reordering A using " << brm->
toString() << std::endl;
114 Teuchos::rcp_dynamic_cast<
const ReorderedBlockedCrsMatrix>(
115 Xpetra::buildReorderedBlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(brm, A));
118 "Block reordering of " << A->Rows() <<
"x" << A->Cols()
119 <<
" blocked operator failed.");
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;
130 Teuchos::rcp_const_cast<ReorderedBlockedCrsMatrix>(brop);
135 currentLevel.Set(
"A", Teuchos::rcp_dynamic_cast<Matrix>(A),
this);
Important warning messages (one line)
Exception indicating invalid cast attempted.
virtual std::string toString() const
Teuchos::RCP< const Xpetra::BlockReorderManager > blockedReorderFromString(std::string reorder)
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Print additional debugging information.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static Teuchos::RCP< Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > SplitMatrix(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &input, Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > rangeMapExtractor, Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > domainMapExtractor, Teuchos::RCP< const Xpetra::MapExtractor< Scalar, LocalOrdinal, GlobalOrdinal, Node > > columnMapExtractor=Teuchos::null, bool bThyraMode=false)
static const NoFactory * get()
RCP< const ParameterList > GetValidParameterList() const
Input.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
static const RCP< const NoFactory > getRCP()
Static Get() functions.
void Build(Level ¤tLevel) const
Build an object with this factory.
void DeclareInput(Level ¤tLevel) const
Specifies the data that this class needs, and the factories that generate that data.
Exception throws to report errors in the internal logical of the program.
Description of what is happening (more verbose)