MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_TransPFactory_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_TRANSPFACTORY_DEF_HPP
11#define MUELU_TRANSPFACTORY_DEF_HPP
12
13#include <Teuchos_ParameterList.hpp>
14#include <Teuchos_Time.hpp>
15
16#include <Xpetra_Matrix.hpp>
17
19
20#include "MueLu_Monitor.hpp"
21#include "MueLu_PerfUtils.hpp"
22#include "MueLu_Utilities.hpp"
23
24namespace MueLu {
25
26template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
28 RCP<ParameterList> validParamList = rcp(new ParameterList());
29 validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Generating factory of the matrix P");
30
31 // Make sure we don't recursively validate options for the matrixmatrix kernels
32 ParameterList norecurse;
33 norecurse.disableRecursiveValidation();
34 validParamList->set<ParameterList>("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
35
36 return validParamList;
37}
38
39template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
41 Input(coarseLevel, "P");
42}
43
44template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
46 FactoryMonitor m(*this, "Transpose P", coarseLevel);
47 std::string label = "MueLu::TransP-" + Teuchos::toString(coarseLevel.GetLevelID());
48
49 RCP<Matrix> P = Get<RCP<Matrix> >(coarseLevel, "P");
50 // If we failed to create a valid P (e.g., # of global aggregates is zero), then we just bail here
51 // This level will ultimately be removed in MueLu_Hierarchy_defs.h via a resize()
52 if (P == Teuchos::null) return;
53
54 const Teuchos::ParameterList& pL = GetParameterList();
55
56 // Reuse pattern if available (multiple solve)
57 RCP<ParameterList> Tparams;
58 if (pL.isSublist("matrixmatrix: kernel params"))
59 Tparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
60 else
61 Tparams = rcp(new ParameterList);
62
63 // By default, we don't need global constants for transpose
64 Tparams->set("compute global constants: temporaries", Tparams->get("compute global constants: temporaries", false));
65 Tparams->set("compute global constants", Tparams->get("compute global constants", false));
66
67 RCP<Matrix> R = Utilities::Transpose(*P, true, label, Tparams);
68
69 if (IsPrint(Statistics2)) {
70 RCP<ParameterList> params = rcp(new ParameterList());
71 params->set("printLoadBalancingInfo", true);
72 params->set("printCommInfo", true);
73 GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*R, "R", params);
74 }
75
76 Set(coarseLevel, "R", R);
77
79 if (P->IsView("stridedMaps"))
80 R->CreateView("stridedMaps", P, true);
82}
83
84} // namespace MueLu
85
86#endif // MUELU_TRANSPFACTORY_DEF_HPP
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.