MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_ShiftedLaplacianOperator_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_SHIFTEDLAPLACIANOPERATOR_DEF_HPP
11#define MUELU_SHIFTEDLAPLACIANOPERATOR_DEF_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14
15#include <Xpetra_Matrix.hpp>
16#include <Xpetra_CrsMatrixWrap.hpp>
17#include <Xpetra_BlockedCrsMatrix.hpp>
18#include <Xpetra_TpetraMultiVector.hpp>
19#include <Xpetra_MultiVectorFactory.hpp>
20
22#include "MueLu_Hierarchy.hpp"
23#include "MueLu_Utilities.hpp"
24
25namespace MueLu {
26
27// ------------- getDomainMap -----------------------
28
29template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> >
32 getDomainMap() const {
33 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XMatrix;
34
35 RCP<MueLu::Level> L0 = Hierarchy_->GetLevel(0);
36 RCP<XMatrix> A = L0->Get<RCP<XMatrix> >("A");
37
38 RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > tpbA =
39 Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(A);
40 if (tpbA != Teuchos::null) {
41 return Xpetra::toTpetraNonZero(tpbA->getDomainMap());
42 }
43
44 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > tpA =
45 toTpetra(A);
46 return tpA->getDomainMap();
47}
48
49// ------------- getRangeMap -----------------------
50
51template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
52Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> >
54 getRangeMap() const {
55 typedef Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> XMatrix;
56
57 RCP<MueLu::Level> L0 = Hierarchy_->GetLevel(0);
58 RCP<XMatrix> A = L0->Get<RCP<XMatrix> >("A");
59
60 RCP<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > tpbA =
61 Teuchos::rcp_dynamic_cast<Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >(A);
62 if (tpbA != Teuchos::null)
63 return Xpetra::toTpetraNonZero(tpbA->getRangeMap());
64
65 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > tpA =
66 toTpetra(A);
67 return tpA->getRangeMap();
68}
69
70// ------------- apply -----------------------
71
72template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
73void ShiftedLaplacianOperator<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X,
74 Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Y,
75 Teuchos::ETransp /* mode */, Scalar /* alpha */, Scalar /* beta */) const {
76 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TMV;
77 typedef Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> XTMV;
78 // typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XMV; // unused
79
80 TMV& temp_x = const_cast<TMV&>(X);
81 const XTMV tX(rcpFromRef(temp_x));
82 XTMV tY(rcpFromRef(Y));
83
84 try {
85 tY.putScalar(0.0);
86 Hierarchy_->Iterate(tX, tY, cycles_, true);
87 }
88
89 catch (std::exception& e) {
90 // FIXME add message and rethrow
91 std::cerr << "Caught an exception in MueLu::ShiftedLaplacianOperator::ApplyInverse():" << std::endl
92 << e.what() << std::endl;
93 }
94
95 // update solution with 2-grid error correction
96 /*if(option_==1) {
97 for(int j=0; j<cycles_; j++) {
98 RCP<XMV> residual = MueLu::Utilities<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Residual(*A_, tY, tX);
99 RCP<XMV> coarseResidual = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(R_->getRangeMap(), tX.getNumVectors());
100 RCP<XMV> coarseError = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(R_->getRangeMap(), tX.getNumVectors());
101 R_ -> apply(*residual, *coarseResidual, Teuchos::NO_TRANS, (Scalar) 1.0, (Scalar) 0.0);
102 RCP<TMV> tcoarseR = MueLu::Utilities<Scalar,LocalOrdinal,GlobalOrdinal,Node>::MV2NonConstTpetraMV(coarseResidual);
103 RCP<TMV> tcoarseE = MueLu::Utilities<Scalar,LocalOrdinal,GlobalOrdinal,Node>::MV2NonConstTpetraMV(coarseError);
104 BelosLP_ -> setProblem(tcoarseE,tcoarseR);
105 BelosSM_ -> solve();
106 RCP<XMV> fineError = Xpetra::MultiVectorFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(P_->getRangeMap(), tX.getNumVectors());
107 XTMV tmpcoarseE(rcpFromRef(*tcoarseE));
108 P_ -> apply(tmpcoarseE, *fineError, Teuchos::NO_TRANS, (Scalar) 1.0, (Scalar) 0.0);
109 tY.update((Scalar) 1.0, *fineError, (Scalar) 1.0);
110 }
111 }
112
113 try {
114 Hierarchy_->Iterate(tX, tY, 1, false);
115 }
116
117 catch(std::exception& e) {
118 //FIXME add message and rethrow
119 std::cerr << "Caught an exception in MueLu::ShiftedLaplacianOperator::ApplyInverse():" << std::endl
120 << e.what() << std::endl;
121 }*/
122}
123
124// ------------- apply -----------------------
125template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
129
130} // namespace MueLu
131
132#endif // ifdef MUELU_SHIFTEDLAPLACIANOPERATOR_DEF_HPP
MueLu::DefaultScalar Scalar
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
void apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::one()) const
Returns in Y the result of a Tpetra::Operator applied to a Tpetra::MultiVector X.
Namespace for MueLu classes and methods.