MueLu  Version of the Day
MueLu_LocalOrdinalTransferFactory_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_LOCALORDINALTRANSFER_FACTORY_DECL_HPP
47 #define MUELU_LOCALORDINALTRANSFER_FACTORY_DECL_HPP
48 
49 #include "MueLu_ConfigDefs.hpp"
51 #include "Xpetra_MultiVector_fwd.hpp"
52 #include "Xpetra_MultiVectorFactory_fwd.hpp"
53 #include "Xpetra_CrsGraph_fwd.hpp"
54 
57 #include "MueLu_Utilities_fwd.hpp"
58 
59 namespace MueLu {
60 
98  template<class LocalOrdinal = DefaultLocalOrdinal,
100  class Node = DefaultNode>
102 #undef MUELU_LOCALORDINALTRANSFERFACTORY_SHORT
104 
105  public:
107 
108 
109  // Default constructor is distabled
110  LocalOrdinalTransferFactory() = delete;
111 
120  LocalOrdinalTransferFactory(const std::string & TransferVecName, const std::string & mode): TransferVecName_(TransferVecName) {
121  if(mode == "classical") useAggregatesMode_ = false;
122  else useAggregatesMode_ = true;
123  }
124 
127 
128  RCP<const ParameterList> GetValidParameterList() const;
129 
131 
133 
134 
140  void DeclareInput(Level &finelevel, Level &coarseLevel) const;
141 
143 
145 
146 
148  void Build(Level & fineLevel, Level &coarseLevel) const;
149 
151 
152  private:
153 
154  void BuildAggregates(Level & fineLevel, Level &coarseLevel) const;
155 
156  void BuildFC(Level & fineLevel, Level &coarseLevel) const;
157 
160 
162  std::string TransferVecName_;
163 
164  }; // class LocalOrdinalTransferFactory
165 
166 } // namespace MueLu
167 
168 #define MUELU_LOCALORDINALTRANSFERFACTORY_SHORT
169 #endif // MUELU_LOCALORDINALTRANSFER_FACTORY_DECL_HPP
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
MueLu::DefaultLocalOrdinal LocalOrdinal
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
bool useAggregatesMode_
Use aggregates mode (as opposed to FC mode)
std::string TransferVecName_
The name for the vector to be transfered. This allows us to have multiple factories for different var...
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Class for transferring a vector of local ordinals from a finer level to a coarser one...
Base class for factories that use two levels (fineLevel and coarseLevel).
Namespace for MueLu classes and methods.
MueLu::DefaultNode Node
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void DeclareInput(Level &finelevel, Level &coarseLevel) const
Specifies the data that this class needs, and the factories that generate that data.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
LocalOrdinalTransferFactory(const std::string &TransferVecName, const std::string &mode)
Constructor.
void BuildFC(Level &fineLevel, Level &coarseLevel) const
void BuildAggregates(Level &fineLevel, Level &coarseLevel) const