MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_PermutingSmoother_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 /*
47  * MueLu_PermutingSmoother_decl.hpp
48  *
49  * Created on: Nov 28, 2012
50  * Author: wiesner
51  */
52 
53 #ifndef MUELU_PERMUTINGSMOOTHER_DECL_HPP
54 #define MUELU_PERMUTINGSMOOTHER_DECL_HPP
55 
57 #include <Xpetra_CrsMatrixWrap_fwd.hpp>
58 #include <Xpetra_Matrix_fwd.hpp>
59 #include <Xpetra_Vector_fwd.hpp>
60 #include <Xpetra_VectorFactory_fwd.hpp>
61 #include <Xpetra_MultiVector_fwd.hpp>
62 #include <Xpetra_MultiVectorFactory_fwd.hpp>
63 
64 #include "MueLu_ConfigDefs.hpp"
65 
66 #include "MueLu_SmootherPrototype.hpp"
69 #include "MueLu_Utilities_fwd.hpp"
71 
72 namespace MueLu {
73 
81  template <class Scalar = SmootherPrototype<>::scalar_type,
82  class LocalOrdinal = typename SmootherPrototype<Scalar>::local_ordinal_type,
83  class GlobalOrdinal = typename SmootherPrototype<Scalar, LocalOrdinal>::global_ordinal_type,
84  class Node = typename SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
85  class PermutingSmoother : public SmootherPrototype<Scalar,LocalOrdinal,GlobalOrdinal,Node>
86  {
87 #undef MUELU_PERMUTINGSMOOTHER_SHORT
88 #include "MueLu_UseShortNames.hpp"
89 
90  public:
91 
93 
94 
102  PermutingSmoother(std::string const & mapName, const RCP<const FactoryBase> & mapFact, std::string const & type = "", const Teuchos::ParameterList& paramList = Teuchos::ParameterList(), LO const & overlap = 0, RCP<FactoryBase> permFact = Teuchos::null);
103 
105  virtual ~PermutingSmoother();
107 
109 
110 
111  void DeclareInput(Level &currentLevel) const;
112 
114 
116 
117 
119  void Setup(Level &currentLevel);
120 
127  void Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero = false) const;
129 
131 
133 
134 
136  std::string description() const;
137 
139  //using MueLu::Describable::describe; // overloading, not hiding
140  void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel = Default) const;
141 
143  size_t getNodeSmootherComplexity() const {return s_->getNodeSmootherComplexity();}
144 
146 
147  private:
149  std::string type_;
150 
153 
156 
159 
162 
165 
166  //
167  // Underlying Smoother
168  //
169 
171  RCP<SmootherPrototype> s_; // TrilinosSmoother object
172 
173  }; // class PermutingSmoother
174 
175 } // namespace MueLu
176 
177 #define MUELU_PERMUTINGSMOOTHER_SHORT
178 #endif /* MUELU_PERMUTINGSMOOTHER_DECL_HPP */
std::string description() const
Return a simple one-line description of this object.
MueLu::DefaultLocalOrdinal LocalOrdinal
RCP< SmootherPrototype > Copy() const
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the direct solver. Solves the linear system AX=B using the constructed solver.
Base class for smoother prototypes.
MueLu::DefaultNode Node
This class first calculates row- and column permutation operators and applies a smoother to the permu...
LO overlap_
overlap when using the smoother in additive Schwarz mode
void Setup(Level &currentLevel)
Set up the direct solver.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Teuchos::RCP< Matrix > diagScalingOp_
scaling matrix object
RCP< FactoryBase > permFact_
Permutation Factory.
void DeclareInput(Level &currentLevel) const
Input.
RCP< Matrix > permP_
permP matrix object
std::string type_
ifpack1/2-specific key phrase that denote smoother type
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
PermutingSmoother(std::string const &mapName, const RCP< const FactoryBase > &mapFact, std::string const &type="", const Teuchos::ParameterList &paramList=Teuchos::ParameterList(), LO const &overlap=0, RCP< FactoryBase > permFact=Teuchos::null)
Constructor.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
RCP< Matrix > permQT_
permQT matrix object
RCP< SmootherPrototype > s_
Smoother.