MueLu  Version of the Day
MueLu_IntrepidPCoarsenFactory_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_IPCFACTORY_DECL_HPP
47 #define MUELU_IPCFACTORY_DECL_HPP
48 
49 #include <string>
50 #include <vector>
51 
52 #include "MueLu_ConfigDefs.hpp"
54 
55 #include "MueLu_Level_fwd.hpp"
57 #include "MueLu_PerfUtils_fwd.hpp"
58 #include "MueLu_PFactory.hpp"
61 #include "MueLu_Utilities_fwd.hpp"
62 
63 #include "Intrepid2_Basis.hpp"
64 
65 #ifdef HAVE_MUELU_INTREPID2_REFACTOR
66 #include "Kokkos_DynRankView.hpp"
67 #else
68 #include "Intrepid2_FieldContainer.hpp"
69 #endif
70 
71 #include <Xpetra_Import.hpp>
72 
73 namespace MueLu {
74 
115  template <class Scalar = double, class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = KokkosClassic::DefaultNode::DefaultNodeType>
117 #undef MUELU_INTREPIDPCOARSENFACTORY_SHORT
118 #include "MueLu_UseShortNames.hpp"
119 
120  public:
121 #ifdef HAVE_MUELU_INTREPID2_REFACTOR
122  typedef Kokkos::DynRankView<LocalOrdinal,typename Node::device_type> LOFieldContainer;
123  typedef Kokkos::DynRankView<double,typename Node::device_type> SCFieldContainer;
124  typedef Intrepid2::Basis<typename Node::device_type::execution_space,double,double> Basis; // Hardwired on purpose
125 #else
126  typedef Intrepid2::FieldContainer<LocalOrdinal> LOFieldContainer;
127  typedef Intrepid2::FieldContainer<double> SCFieldContainer;
128  typedef Intrepid2::Basis<double,SCFieldContainer>
129 #endif
130 
132 
137 
140 
141  RCP<const ParameterList> GetValidParameterList() const;
142 
144 
146 
147 
148  void DeclareInput(Level &fineLevel, Level &coarseLevel) const;
149 
151 
153 
154 
160  void Build(Level& fineLevel, Level &coarseLevel) const;
161 
162  void BuildP(Level &fineLevel, Level &coarseLevel) const; //Build()
163 
165  private:
167 
168  // NOTE: This is hardwired to double on purpose.
169  void GenerateLinearCoarsening_pn_kirby_to_p1(const LOFieldContainer & hi_elemToNode,
170  const std::vector<bool> & hi_nodeIsOwned,
171  const SCFieldContainer & hi_DofCoords,
172  const std::vector<size_t> &lo_node_in_hi,
173  const Basis &lo_Basis,
174  const std::vector<LocalOrdinal> & hi_to_lo_map,
175  const Teuchos::RCP<const Map> & lo_colMap,
176  const Teuchos::RCP<const Map> & lo_domainMap,
177  const Teuchos::RCP<const Map> & hi_map,
178  Teuchos::RCP<Matrix>& P) const;
180 
181  // NOTE: This is hardwired to double on purpose.
182  void GenerateLinearCoarsening_pn_kirby_to_pm(const LOFieldContainer & hi_elemToNode,
183  const std::vector<bool> & hi_nodeIsOwned,
184  const SCFieldContainer & hi_DofCoords,
185  const LOFieldContainer & lo_elemToHiRepresentativeNode,
186  const Basis &lo_basis,
187  const std::vector<LocalOrdinal> & hi_to_lo_map,
188  const Teuchos::RCP<const Map> & lo_colMap,
189  const Teuchos::RCP<const Map> & lo_domainMap,
190  const Teuchos::RCP<const Map> & hi_map,
191  Teuchos::RCP<Matrix>& P) const;
192 
193 
194  }; //class IntrepidPCoarsenFactory
195 
196 
197  /* Utility functions for use with Intrepid */
198  namespace MueLuIntrepid {
199 
200 #ifdef HAVE_MUELU_INTREPID2_REFACTOR
201  template<class Scalar,class KokkosExecutionSpace>
202  Teuchos::RCP<Intrepid2::Basis<KokkosExecutionSpace,Scalar,Scalar> > BasisFactory(const std::string & name, int & degree);
203 #else
204  // NOTE: This function will not work with Stokhos scalar types, due to deficiencies upstream.
205  template<class Scalar>
206  Teuchos::RCP<Intrepid2::Basis<Scalar,Intrepid2::FieldContainer<Scalar> > > BasisFactory(const std::string & name, int & degree);
207 #endif
208 
209 #ifdef HAVE_MUELU_INTREPID2_REFACTOR
210  template<class Scalar,class KokkosDeviceType>
211  void IntrepidGetLoNodeInHi(const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space,Scalar,Scalar> > &hi_basis,
212  const Teuchos::RCP<Intrepid2::Basis<typename KokkosDeviceType::execution_space,Scalar,Scalar> > &lo_basis,
213  std::vector<size_t> & lo_node_in_hi,
214  Kokkos::DynRankView<Scalar,KokkosDeviceType> & hi_DofCoords);
215 #else
216  template <class Scalar, class ArrayScalar>
217  void IntrepidGetLoNodeInHi(const Teuchos::RCP<Intrepid2::Basis<Scalar,ArrayScalar> > &hi_basis,
218  const Teuchos::RCP<Intrepid2::Basis<Scalar,ArrayScalar> > &lo_basis,
219  std::vector<size_t> & lo_node_in_hi,
220  ArrayScalar & hi_DofCoords);
221 #endif
222 
223  template<class LocalOrdinal, class GlobalOrdinal, class Node, class LOFieldContainer>
224  void GenerateLoNodeInHiViaGIDs(const std::vector<std::vector<size_t> > & candidates,const LOFieldContainer & hi_elemToNode,
225  RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & hi_columnMap,
226  LOFieldContainer & lo_elemToHiRepresentativeNode);
227 
228  template <class LocalOrdinal, class LOFieldContainer>
229  void BuildLoElemToNode(const LOFieldContainer & hi_elemToNode,
230  const std::vector<bool> & hi_nodeIsOwned,
231  const std::vector<size_t> & lo_node_in_hi,
232  const Teuchos::ArrayRCP<const int> & hi_isDirichlet,
233  LOFieldContainer & lo_elemToNode,
234  std::vector<bool> & lo_nodeIsOwned,
235  std::vector<LocalOrdinal> & hi_to_lo_map,
236  int & lo_numOwnedNodes);
237 
238  template <class LocalOrdinal, class LOFieldContainer>
239  void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer & hi_elemToNode,
240  const std::vector<bool> & hi_nodeIsOwned,
241  const LOFieldContainer & lo_elemToHiRepresentativeNode,
242  LOFieldContainer & lo_elemToNode,
243  std::vector<bool> & lo_nodeIsOwned,
244  std::vector<LocalOrdinal> & hi_to_lo_map,
245  int & lo_numOwnedNodes);
246 
247 
248  template <class LocalOrdinal, class GlobalOrdinal, class Node>
249  void GenerateColMapFromImport(const Xpetra::Import<LocalOrdinal,GlobalOrdinal,Node> & hi_importer,const std::vector<LocalOrdinal> &hi_to_lo_map,const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> & lo_domainMap, const size_t & lo_columnMapLength, RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > & lo_columnMap);
250 
251 
252  template<class Basis, class SCFieldContainer>
253  void GenerateRepresentativeBasisNodes(const Basis & basis, const SCFieldContainer & ReferenceNodeLocations, const double threshold, std::vector<std::vector<size_t> > & representative_node_candidates);
254 
255  // ! Given an element to (global) node map and a basis, determine one global ordinal per geometric entity (vertex, edge, face,
256  // ! interior). On exit, seeds container is of dimension (spaceDim+1), and contains a sorted vector of local ordinals
257  // ! belonging to entities of that dimension. Only locally-owned degrees of freedom (as determined by rowMap and columnMap)
258  // ! will be stored in seeds.
259  template<class Basis, class LOFieldContainer, class LocalOrdinal, class GlobalOrdinal, class Node>
260  void FindGeometricSeedOrdinals(Teuchos::RCP<Basis> basis, const LOFieldContainer &elementToNodeMap,
261  std::vector<std::vector<LocalOrdinal> > &seeds,
262  const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &rowMap,
263  const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> &columnMap);
264 
265  }//namespace MueLuIntrepid
266 } //namespace MueLu
267 
268 #define MUELU_INTREPIDPCOARSENFACTORY_SHORT
269 #endif // MUELU_IPCFACTORY_DECL_HPP
void GenerateLoNodeInHiViaGIDs(const std::vector< std::vector< size_t > > &candidates, const LOFieldContainer &hi_elemToNode, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &hi_columnMap, LOFieldContainer &lo_elemToHiRepresentativeNode)
void FindGeometricSeedOrdinals(Teuchos::RCP< Basis > basis, const LOFieldContainer &elementToNodeMap, std::vector< std::vector< LocalOrdinal > > &seeds, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &rowMap, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &columnMap)
void IntrepidGetLoNodeInHi(const Teuchos::RCP< Intrepid2::Basis< Scalar, ArrayScalar > > &hi_basis, const Teuchos::RCP< Intrepid2::Basis< Scalar, ArrayScalar > > &lo_basis, std::vector< size_t > &lo_node_in_hi, ArrayScalar &hi_DofCoords)
Factory for building transfer operators based on coarsening in polynomial degree, following the Intre...
void GenerateLinearCoarsening_pn_kirby_to_pm(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const LOFieldContainer &lo_elemToHiRepresentativeNode, const Basis &lo_basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
Namespace for MueLu classes and methods.
Intrepid2::FieldContainer< double > SCFieldContainer
IntrepidPCoarsenFactory()
Constructor. User can supply a factory for generating the tentative prolongator.
Teuchos::RCP< Intrepid2::Basis< Scalar, Intrepid2::FieldContainer< Scalar > > > BasisFactory(const std::string &name, int &degree)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void BuildLoElemToNode(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const std::vector< size_t > &lo_node_in_hi, const Teuchos::ArrayRCP< const int > &hi_isDirichlet, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
void BuildLoElemToNodeViaRepresentatives(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const LOFieldContainer &lo_elemToHiRepresentativeNode, LOFieldContainer &lo_elemToNode, std::vector< bool > &lo_nodeIsOwned, std::vector< LocalOrdinal > &hi_to_lo_map, int &lo_numOwnedNodes)
void GenerateLinearCoarsening_pn_kirby_to_p1(const LOFieldContainer &hi_elemToNode, const std::vector< bool > &hi_nodeIsOwned, const SCFieldContainer &hi_DofCoords, const std::vector< size_t > &lo_node_in_hi, const Basis &lo_Basis, const std::vector< LocalOrdinal > &hi_to_lo_map, const Teuchos::RCP< const Map > &lo_colMap, const Teuchos::RCP< const Map > &lo_domainMap, const Teuchos::RCP< const Map > &hi_map, Teuchos::RCP< Matrix > &P) const
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Intrepid2::FieldContainer< LocalOrdinal > LOFieldContainer
void GenerateColMapFromImport(const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &hi_importer, const std::vector< LocalOrdinal > &hi_to_lo_map, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &lo_domainMap, const size_t &lo_columnMapLength, RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &lo_columnMap)
void GenerateRepresentativeBasisNodes(const Basis &basis, const SCFieldContainer &ReferenceNodeLocations, const double threshold, std::vector< std::vector< size_t > > &representative_node_candidates)
Factory that provides an interface for a concrete implementation of a prolongation operator...
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.