10#ifndef MUELU_ITERATOROPS_HPP_
11#define MUELU_ITERATOROPS_HPP_
15#include "Xpetra_Matrix.hpp"
16#include "Xpetra_MatrixMatrix.hpp"
17#include "Xpetra_CrsMatrixWrap.hpp"
19#include "Xpetra_Map.hpp"
20#include "Xpetra_StridedMap.hpp"
21#include "Xpetra_StridedMapFactory.hpp"
22#include "Xpetra_MapExtractor.hpp"
23#include "Xpetra_MatrixFactory.hpp"
24#include "Xpetra_BlockedCrsMatrix.hpp"
32template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
35 const Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Dinv,
36 const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
37 const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& B,
38 Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& C,
39 bool call_FillComplete_on_result =
true,
40 bool doOptimizeStorage =
true,
41 const std::string& label = std::string(),
42 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
49 "MueLu::Jacobi: row map of C is not same as row map of A")
51 "MueLu::Jacobi: row map of C is not same as row map of B");
55 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
57 if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
58#ifndef HAVE_MUELU_EPETRAEXT
63 }
else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
64 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpA = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(A);
65 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpB = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(B);
66 Tpetra::CrsMatrix<SC, LO, GO, NO>& tpC = Xpetra::Helpers<SC, LO, GO, NO>::Op2NonConstTpetraCrs(C);
67 const RCP<Tpetra::Vector<SC, LO, GO, NO> >& tpD = toTpetra(Dinv);
68 Tpetra::MatrixMatrix::Jacobi(omega, *tpD, tpA, tpB, tpC, haveMultiplyDoFillComplete, label, params);
71 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
72 RCP<Teuchos::ParameterList> fillParams = rcp(
new Teuchos::ParameterList());
73 fillParams->set(
"Optimize Storage", doOptimizeStorage);
74 C.fillComplete(B.getDomainMap(), B.getRangeMap(), fillParams);
78 RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpA = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(A));
79 RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpB = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(B));
80 C.CreateView(
"stridedMaps", rcpA,
false, rcpB,
false);
83#if defined(HAVE_MUELU_EPETRA) && !defined(XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES)
86 const Xpetra::Vector<double, int, int, Xpetra::EpetraNode>& Dinv,
87 const Xpetra::Matrix<double, int, int, Xpetra::EpetraNode>& A,
88 const Xpetra::Matrix<double, int, int, Xpetra::EpetraNode>& B,
89 Xpetra::Matrix<double, int, int, Xpetra::EpetraNode>& C,
90 bool call_FillComplete_on_result,
91 bool doOptimizeStorage,
92 const std::string& label,
93 const Teuchos::RCP<Teuchos::ParameterList>& params);
96#if defined(HAVE_MUELU_EPETRA) && !defined(XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES)
99 const Xpetra::Vector<double, int, long long, Xpetra::EpetraNode>& Dinv,
100 const Xpetra::Matrix<double, int, long long, Xpetra::EpetraNode>& A,
101 const Xpetra::Matrix<double, int, long long, Xpetra::EpetraNode>& B,
102 Xpetra::Matrix<double, int, long long, Xpetra::EpetraNode>& C,
103 bool call_FillComplete_on_result,
104 bool doOptimizeStorage,
105 const std::string& label,
106 const Teuchos::RCP<Teuchos::ParameterList>& params);
116template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
119 static RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >
120 Jacobi(
Scalar omega,
const Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Dinv,
const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& A,
const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& B, RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > C_in, Teuchos::FancyOStream& fos,
const std::string& label, RCP<ParameterList>& params) {
124 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > C = C_in;
125 if (C == Teuchos::null) {
126 C = Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(B.getRowMap(), Teuchos::OrdinalTraits<LocalOrdinal>::zero());
129 fos <<
"Reuse C pattern" << std::endl;
132 MueLu::Jacobi<Scalar, LocalOrdinal, GlobalOrdinal, Node>(omega, Dinv, A, B, *C,
true,
true, label, params);
133 C->CreateView(
"stridedMaps", rcpFromRef(A),
false, rcpFromRef(B),
false);
141#define MUELU_ITERATOROPS_SHORT
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception throws to report errors in the internal logical of the program.
MueLu utility class containing iteration operators.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Jacobi(Scalar omega, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Dinv, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > C_in, Teuchos::FancyOStream &fos, const std::string &label, RCP< ParameterList > ¶ms)
Namespace for MueLu classes and methods.
void Jacobi< double, int, long long, Xpetra::EpetraNode >(double omega, const Xpetra::Vector< double, int, long long, Xpetra::EpetraNode > &Dinv, const Xpetra::Matrix< double, int, long long, Xpetra::EpetraNode > &A, const Xpetra::Matrix< double, int, long long, Xpetra::EpetraNode > &B, Xpetra::Matrix< double, int, long long, Xpetra::EpetraNode > &C, bool call_FillComplete_on_result, bool doOptimizeStorage, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
void Jacobi(Scalar omega, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Dinv, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &C, bool call_FillComplete_on_result=true, bool doOptimizeStorage=true, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
void Jacobi< double, int, int, Xpetra::EpetraNode >(double omega, const Xpetra::Vector< double, int, int, Xpetra::EpetraNode > &Dinv, const Xpetra::Matrix< double, int, int, Xpetra::EpetraNode > &A, const Xpetra::Matrix< double, int, int, Xpetra::EpetraNode > &B, Xpetra::Matrix< double, int, int, Xpetra::EpetraNode > &C, bool call_FillComplete_on_result, bool doOptimizeStorage, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)