10 #ifndef MUELU_MULTIPHYS_DEF_HPP
11 #define MUELU_MULTIPHYS_DEF_HPP
16 #include "Xpetra_Map.hpp"
21 #include "MueLu_SaPFactory.hpp"
22 #include "MueLu_AggregationExportFactory.hpp"
23 #include "MueLu_Utilities.hpp"
25 #include "MueLu_Hierarchy.hpp"
26 #include "MueLu_RAPFactory.hpp"
28 #include "MueLu_PerfUtils.hpp"
29 #include "MueLu_ParameterListInterpreter.hpp"
31 #include <MueLu_HierarchyUtils.hpp>
36 #ifdef HAVE_MUELU_CUDA
37 #include "cuda_profiler_api.h"
42 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
44 return AmatMultiphysics_->getDomainMap();
47 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
49 return AmatMultiphysics_->getRangeMap();
52 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
56 list.
set(
"multigrid algorithm",
"combine");
57 list.
set(
"combine: numBlks", nBlks_);
60 std::string verbosity = list.
get(
"verbosity",
"high");
63 arrayOfParamLists_.resize(nBlks_);
64 for (
int i = 0; i < nBlks_; i++) {
67 arrayOfParamLists_[i] = Teuchos::rcpFromRef(list.
sublist(listName));
71 arrayOfParamLists_[i]->set(
"verbosity", arrayOfParamLists_[i]->
get(
"verbosity", verbosity));
72 arrayOfParamLists_[i]->set(
"smoother: pre or post",
"none");
73 arrayOfParamLists_[i]->set(
"smoother: type",
"none");
77 useKokkos_ = !Node::is_serial;
78 useKokkos_ = list.
get(
"use kokkos refactor", useKokkos_);
80 paramListMultiphysics_ = Teuchos::rcpFromRef(list);
83 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
96 std::string timerLabel;
98 timerLabel =
"MueLu MultiPhys: compute (reuse)";
100 timerLabel =
"MueLu MultiPhys: compute";
106 for (
int iii = 0; iii < nBlks_; iii++) {
107 if (arrayOfCoords_ != Teuchos::null) {
108 arrayOfParamLists_[iii]->sublist(
"user data").set(
"Coordinates", arrayOfCoords_[iii]);
111 bool wantToRepartition =
false;
112 if (paramListMultiphysics_->isParameter(
"repartition: enable"))
113 wantToRepartition = paramListMultiphysics_->get<
bool>(
"repartition: enable");
115 arrayOfParamLists_[iii]->set(
"repartition: enable", wantToRepartition);
116 arrayOfParamLists_[iii]->set(
"repartition: rebalance P and R",
true);
117 arrayOfParamLists_[iii]->set(
"repartition: explicit via new copy rebalance P and R",
true);
119 if (paramListMultiphysics_->isParameter(
"repartition: use subcommunicators"))
120 arrayOfParamLists_[iii]->set(
"repartition: use subcommunicators", paramListMultiphysics_->isParameter(
"repartition: use subcommunicators"));
122 arrayOfParamLists_[iii]->set(
"repartition: use subcommunicators",
true);
127 paramListMultiphysics_->set<
bool>(
"repartition: enable",
false);
130 for (
int i = 0; i < nBlks_; i++) {
132 arrayOfAuxMatrices_[i]->setObjectLabel(operatorLabel);
134 LO tempNlevels = arrayOfHierarchies_[i]->GetGlobalNumLevels();
135 if (tempNlevels < maxLevels) maxLevels = tempNlevels;
139 for (
LO i = 0; i < maxLevels; i++) {
140 hierarchyMultiphysics_->AddNewLevel();
142 for (
int i = 0; i < nBlks_; i++) {
146 paramListMultiphysics_->set(
"coarse: max size", 1);
147 paramListMultiphysics_->set(
"max levels", maxLevels);
159 const std::string& levelName = inListEntry->first;
160 if (levelName.find(
"subblockList") != 0) stripped.
setEntry(inListEntry->first, inListEntry->second);
164 hierarchyMultiphysics_->setlib(AmatMultiphysics_->getDomainMap()->lib());
165 hierarchyMultiphysics_->SetProcRankVerbose(AmatMultiphysics_->getDomainMap()->getComm()->getRank());
168 hierarchyMultiphysics_->GetLevel(0)->Set(
"A", AmatMultiphysics_);
173 mueLuFactory->SetupHierarchy(*hierarchyMultiphysics_);
178 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
190 return Teuchos::null;
193 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
195 bool reuse = !AmatMultiphysics_.is_null();
196 AmatMultiphysics_ = AmatMultiphysics_new;
197 if (ComputePrec) compute(reuse);
200 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
202 hierarchyMultiphysics_->Iterate(RHS, X, 1,
true);
205 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
211 hierarchyMultiphysics_->Iterate(RHS, X, 1,
true);
214 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
219 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
227 arrayOfHierarchies_.resize(nBlks_);
228 for (
int i = 0; i < nBlks_; i++) arrayOfHierarchies_[i] = Teuchos::null;
232 enable_reuse_ =
false;
240 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
243 std::ostringstream oss;
247 oss <<
"\n--------------------------------------------------------------------------------\n"
248 <<
"--- MultiPhysics Summary ---\n"
249 "--------------------------------------------------------------------------------"
256 AmatMultiphysics_->getRowMap()->getComm()->barrier();
258 for (
int i = 0; i < nBlks_; i++) {
259 numRows = arrayOfAuxMatrices_[i]->getGlobalNumRows();
260 nnz = arrayOfAuxMatrices_[i]->getGlobalNumEntries();
274 oss <<
"block " << std::setw(rowspacer) <<
" rows " << std::setw(nnzspacer) <<
" nnz " << std::setw(9) <<
" nnz/row" << std::endl;
275 oss <<
"(" <<
Teuchos::toString(i) <<
", " <<
Teuchos::toString(i) <<
")" << std::setw(rowspacer) << numRows << std::setw(nnzspacer) << nnz << std::setw(9) << as<double>(nnz) / numRows << std::endl;
281 for (
int i = 0; i < nBlks_; i++) {
282 arrayOfHierarchies_[i]->describe(out, GetVerbLevel());
289 #define MUELU_MULTIPHYS_SHORT
290 #endif // ifdef MUELU_MULTIPHYS_DEF_HPP
ParameterList & setEntry(const std::string &name, U &&entry)
Various adapters that will create a MueLu preconditioner that is an Xpetra::Matrix.
ConstIterator end() const
static void AddNonSerializableDataToHierarchy(HierarchyManager &HM, Hierarchy &H, const ParameterList &nonSerialList)
Add non-serializable data to Hierarchy.
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void setParameters(Teuchos::ParameterList &list)
Set parameters.
MsgType toVerbLevel(const std::string &verbLevelStr)
static RCP< Time > getNewTimer(const std::string &name)
Teuchos::RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > CreateXpetraPreconditioner(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> op, const Teuchos::ParameterList &inParamList)
Helper function to create a MueLu preconditioner that can be used by Xpetra.Given an Xpetra::Matrix...
static void SetDefaultVerbLevel(const VerbLevel defaultVerbLevel)
Set the default (global) verbosity level.
void applyInverse(const MultiVector &RHS, MultiVector &X) const
apply standard MultiPhys cycle
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const Teuchos::RCP< const Map > getDomainMap() const
Returns the Xpetra::Map object associated with the domain of this operator.
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
const Teuchos::RCP< const Map > getRangeMap() const
Returns the Xpetra::Map object associated with the range of this operator.
bool isSublist(const std::string &name) const
params_t::ConstIterator ConstIterator
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const
ConstIterator begin() const
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< Teuchos::TimeMonitor > getTimer(std::string name, RCP< const Teuchos::Comm< int >> comm=Teuchos::null) const
get a (synced) timer
Print all timing information.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Exception throws to report errors in the internal logical of the program.
void initialize(const Teuchos::RCP< Matrix > &AmatMultiPhysics, const Teuchos::ArrayRCP< RCP< Matrix >> arrayOfAuxMatrices, const Teuchos::ArrayRCP< Teuchos::RCP< MultiVector >> arrayOfNullspaces, const Teuchos::ArrayRCP< Teuchos::RCP< RealValuedMultiVector >> arrayOfCoords, const int nBlks, Teuchos::ParameterList &List)
void resetMatrix(Teuchos::RCP< Matrix > SM_Matrix_new, bool ComputePrec=true)
Reset system matrix.
void compute(bool reuse=false)
Setup the preconditioner.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
std::string toString(const T &t)
long ExtractNonSerializableData(const Teuchos::ParameterList &inList, Teuchos::ParameterList &serialList, Teuchos::ParameterList &nonSerialList)
Extract non-serializable data from level-specific sublists and move it to a separate parameter list...
static void CopyBetweenHierarchies(Hierarchy &fromHierarchy, Hierarchy &toHierarchy, const std::string fromLabel, const std::string toLabel, const std::string dataType)