MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_RAPShiftFactory_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_RAPSHIFTFACTORY_DECL_HPP
47 #define MUELU_RAPSHIFTFACTORY_DECL_HPP
48 
49 #include <string>
50 
51 #include <Xpetra_Matrix_fwd.hpp>
53 #include <Xpetra_Vector_fwd.hpp>
55 
56 #include "MueLu_ConfigDefs.hpp"
58 
60 #include "MueLu_Level_fwd.hpp"
61 #include "MueLu_PerfUtils_fwd.hpp"
63 
64 namespace MueLu {
72 template <class Scalar = DefaultScalar,
75  class Node = DefaultNode>
77 #undef MUELU_RAPSHIFTFACTORY_SHORT
78 #include "MueLu_UseShortNames.hpp"
79 
80  public:
82 
83 
85 
86  virtual ~RAPShiftFactory() {}
87 
89 
91 
92 
94 
95  void DeclareInput(Level &fineLevel, Level &coarseLevel) const;
96 
98 
100 
101  void Build(Level &fineLevel, Level &coarseLevel) const;
103 
105 
106 
108  void SetImplicitTranspose(bool const &implicit) {
109  implicitTranspose_ = implicit;
110  }
111 
112  void SetShifts(std::vector<Scalar> &shifts) {
113  shifts_.clear();
114  shifts_ = shifts;
115  }
116 
118 
120 
125  void AddTransferFactory(const RCP<const FactoryBase> &factory);
126 
127  // TODO add a function to remove a specific transfer factory?
128 
130  size_t NumTransferFactories() const { return transferFacts_.size(); }
131 
133 
134  private:
137 
139  std::vector<RCP<const FactoryBase> > transferFacts_;
140 
141  // vector of shifting terms
142  std::vector<Scalar> shifts_;
143 
144 }; // class RAPShiftFactory
145 
146 } // namespace MueLu
147 
148 #define MUELU_RAPSHIFTFACTORY_SHORT
149 #endif // MUELU_RAPSHIFTFACTORY_DECL_HPP
size_t NumTransferFactories() const
Returns number of transfer factories.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
MueLu::DefaultLocalOrdinal LocalOrdinal
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Factory for building coarse grid matrices, when the matrix is of the form K+a*M. Useful when you want...
Base class for factories that use two levels (fineLevel and coarseLevel).
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
bool implicitTranspose_
If true, the action of the restriction operator action is implicitly defined by the transpose of the ...
MueLu::DefaultNode Node
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
std::vector< RCP< const FactoryBase > > transferFacts_
list of user-defined transfer Factories
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
void SetShifts(std::vector< Scalar > &shifts)
void SetImplicitTranspose(bool const &implicit)
Indicate that the restriction operator action should be implicitly defined by the transpose of the pr...