MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_ShiftedLaplacianOperator_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_SHIFTEDLAPLACIANOPERATOR_DECL_HPP
47 #define MUELU_SHIFTEDLAPLACIANOPERATOR_DECL_HPP
48 
49 #include "MueLu_ConfigDefs.hpp"
50 
51 #include <Tpetra_Operator.hpp>
53 #include "MueLu_Level.hpp"
54 #include "MueLu_Hierarchy_decl.hpp"
55 #include "MueLu_Utilities.hpp"
56 
57 // TODO: Kokkos headers
58 
59 namespace MueLu {
60 
64 template <class Scalar = Tpetra::Operator<>::scalar_type,
69  : public Tpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
75 
76  public:
78 
79 
82  : Hierarchy_(H)
83  , option_(0) {}
84 
87  const RCP<Matrix> A, int cycles, int iters, int option, double tol)
88  : Hierarchy_(H)
89  , A_(A)
90  , cycles_(cycles)
91  , iters_(iters)
92  , option_(option)
93  , tol_(tol) {
94  // setup 2-level correction
95  /*RCP< MueLu::Level > Level1 = H -> GetLevel(1);
96  R_ = Level1 -> Get< RCP<Matrix> >("R");
97  P_ = Level1 -> Get< RCP<Matrix> >("P");
98  //RCP<Matrix> AP = Level1 -> Get< RCP<Matrix> >("AP graph");
99  RCP<Matrix> AP;
100  AP = MUtils::Multiply(*A_, false, *P_, false, AP);
101  // Optimization storage option. If matrix is not changing later, allow this.
102  bool doOptimizedStorage = true;
103  // Reuse coarse matrix memory if available (multiple solve)
104  //RCP<Matrix> Ac = Level1 -> Get< RCP<Matrix> >("RAP graph");
105  RCP<Matrix> Ac;
106  Ac = MUtils::Multiply(*R_, false, *AP, false, Ac, true, doOptimizedStorage);
107  Ac_ = MUtils::Op2NonConstTpetraCrs(Ac);
108 
109  // Setup Belos for two-level correction
110  BelosList_ = rcp( new Teuchos::ParameterList("GMRES") );
111  BelosList_ -> set("Maximum Iterations", iters_ );
112  BelosList_ -> set("Convergence Tolerance", tol_ );
113  BelosLP_ = rcp( new Belos::LinearProblem<Scalar,MV,OP> );
114  BelosLP_ -> setOperator ( Ac_ );
115  BelosSM_ = rcp( new Belos::BlockGmresSolMgr<Scalar,MV,OP>(BelosLP_, BelosList_) );*/
116  }
117 
120 
122 
125 
128 
130 
140 
142  bool hasTransposeApply() const;
143 
144  private:
149 
150  // RCP< Belos::LinearProblem<Scalar,MV,OP> > BelosLP_;
151  // RCP< Belos::SolverManager<Scalar,MV,OP> > BelosSM_;
152 
153  // cycles -> number of V-cycles
154  // iters -> number of GMRES iterations per correction
155  // option -> 0 if no correction is desired
157  double tol_;
158 };
159 
160 } // namespace MueLu
161 
162 #endif // MUELU_SHIFTEDLAPLACIANOPERATOR_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
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. ...
MueLu::DefaultNode Node
Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > OP
bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > R_
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > MV
Wraps an existing MueLu::Hierarchy as a Tpetra::Operator, with an optional two-level correction...
RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Hierarchy_
MueLu::Utilities< Scalar, LocalOrdinal, GlobalOrdinal, Node > MUtils
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > P_
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A_
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Ac_
ShiftedLaplacianOperator(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H)
Constructor.
ShiftedLaplacianOperator(const RCP< MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &H, const RCP< Matrix > A, int cycles, int iters, int option, double tol)
Auxiliary Constructor.
MueLu utility class.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.