MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_TpetraOperator_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_TPETRAOPERATOR_DEF_HPP
11 #define MUELU_TPETRAOPERATOR_DEF_HPP
12 
13 #include "MueLu_ConfigDefs.hpp"
14 
15 #include <Xpetra_BlockedMap.hpp>
16 #include <Xpetra_Matrix.hpp>
17 #include <Xpetra_CrsMatrixWrap.hpp>
19 #include <Xpetra_Operator.hpp>
20 #include <Xpetra_TpetraMultiVector.hpp>
21 
23 #include "MueLu_Hierarchy.hpp"
24 #include "MueLu_Utilities.hpp"
25 
26 namespace MueLu {
27 
28 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
30 
31 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37 
38  RCP<const Map> domainMap;
39  if (!Hierarchy_.is_null())
40  domainMap = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >("A")->getDomainMap();
41  else
42  domainMap = Operator_->getDomainMap();
43 
44  RCP<const BlockedMap> bDomainMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(domainMap);
45  if (bDomainMap.is_null() == false) {
46  return Xpetra::toTpetraNonZero(bDomainMap->getFullMap());
47  }
48  return Xpetra::toTpetraNonZero(domainMap);
49 }
50 
51 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
56 
57  RCP<const Map> rangeMap;
58  if (!Hierarchy_.is_null())
59  rangeMap = Hierarchy_->GetLevel(0)->template Get<RCP<Matrix> >("A")->getRangeMap();
60  else
61  rangeMap = Operator_->getRangeMap();
62 
63  RCP<const BlockedMap> bRangeMap = Teuchos::rcp_dynamic_cast<const BlockedMap>(rangeMap);
64  if (bRangeMap.is_null() == false) {
65  return Xpetra::toTpetraNonZero(bRangeMap->getFullMap());
66  }
67  return Xpetra::toTpetraNonZero(rangeMap);
68 }
69 
70 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
73  Teuchos::ETransp mode, Scalar /* alpha */, Scalar /* beta */) const {
76 
77  TEUCHOS_TEST_FOR_EXCEPTION(mode != Teuchos::NO_TRANS, std::logic_error, "MueLu::TpetraOperator does not support applying the adjoint operator");
78 
79  try {
80  TMV& temp_x = const_cast<TMV&>(X);
81  const XTMV tX(rcpFromRef(temp_x));
82  XTMV tY(rcpFromRef(Y));
83 
84  if (!Hierarchy_.is_null())
85  Hierarchy_->Iterate(tX, tY, 1, true);
86  else
87  Operator_->apply(tX, tY);
88 
89  } catch (std::exception& e) {
90  std::cerr << "MueLu::TpetraOperator::apply : detected an exception" << std::endl
91  << e.what() << std::endl;
92  throw;
93  }
94 }
95 
96 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98  return false;
99 }
100 
101 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
104  return Hierarchy_;
105 }
106 
107 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
110  return Operator_;
111 }
112 
113 } // namespace MueLu
114 
115 #endif // ifdef MUELU_TPETRAOPERATOR_DEF_HPP
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > toTpetraNonZero(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
virtual ~TpetraOperator()
Destructor.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
MueLu::DefaultScalar Scalar
RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetOperator() const
Direct access to the underlying MueLu::Operator.
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetHierarchy() const
Direct access to the underlying MueLu::Hierarchy.
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. ...
bool is_null() const