Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_IntrepidFieldPattern.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef __Panzer_IntrepidFieldPattern_hpp__
12 #define __Panzer_IntrepidFieldPattern_hpp__
13 
14 #include "Panzer_FieldPattern.hpp"
15 
16 // Trilinos includes
17 #include "Kokkos_Core.hpp"
18 #include "Kokkos_DynRankView.hpp"
19 #include "Intrepid2_Basis.hpp"
20 
21 #include "Phalanx_KokkosDeviceTypes.hpp"
22 #include "Teuchos_RCP.hpp"
23 #include <set>
24 
25 namespace panzer {
26 
31  public:
32  Intrepid2FieldPattern(const Teuchos::RCP< Intrepid2::Basis<PHX::Device,double,double> > &intrepidBasis);
33 
34  virtual int getSubcellCount(int dim) const;
35  virtual const std::vector<int> & getSubcellIndices(int dim, int cellIndex) const;
36  virtual int getDimension() const;
37  virtual shards::CellTopology getCellTopology() const;
38 
39  virtual void getSubcellClosureIndices(int dim,int cellIndex,std::vector<int> & indices) const;
40 
41  // static functions for examining shards objects
42 
54  static void buildSubcellClosure(const shards::CellTopology & cellTopo,unsigned dim,unsigned subCell,
55  std::set<std::pair<unsigned,unsigned> > & closure);
56 
68  static void findContainedSubcells(const shards::CellTopology & cellTopo,unsigned dim,
69  const std::vector<unsigned> & nodes,
70  std::set<std::pair<unsigned,unsigned> > & subCells);
71 
79  static void getSubcellNodes(const shards::CellTopology & cellTopo,unsigned dim,unsigned subCell,
80  std::vector<unsigned> & nodes);
81 
90 
96  void getInterpolatoryCoordinates(Kokkos::DynRankView<double,PHX::Device> & coords) const;
97 
107  void getInterpolatoryCoordinates(const Kokkos::DynRankView<double,PHX::Device> & cellNodes,
108  Kokkos::DynRankView<double,PHX::Device> & coords,
109  Teuchos::RCP<const shards::CellTopology> meshCellTopology=Teuchos::null) const;
110 
113 
114  protected:
116 
117  //mutable std::vector<int> subcellIndices_;
118  mutable std::vector<std::vector<std::vector<int> > > subcellIndicies_;
119  std::vector<int> empty_;
120  };
121 
122 }
123 
124 #endif
virtual const std::vector< int > & getSubcellIndices(int dim, int cellIndex) const
static void buildSubcellClosure(const shards::CellTopology &cellTopo, unsigned dim, unsigned subCell, std::set< std::pair< unsigned, unsigned > > &closure)
virtual shards::CellTopology getCellTopology() const
bool supportsInterpolatoryCoordinates() const
Does this field pattern support interpolatory coordinates?
Intrepid2FieldPattern(const Teuchos::RCP< Intrepid2::Basis< PHX::Device, double, double > > &intrepidBasis)
static void getSubcellNodes(const shards::CellTopology &cellTopo, unsigned dim, unsigned subCell, std::vector< unsigned > &nodes)
virtual void getSubcellClosureIndices(int dim, int cellIndex, std::vector< int > &indices) const
Teuchos::RCP< Intrepid2::Basis< PHX::Device, double, double > > getIntrepidBasis() const
Returns the underlying Intrepid2::Basis object.
std::vector< std::vector< std::vector< int > > > subcellIndicies_
static void findContainedSubcells(const shards::CellTopology &cellTopo, unsigned dim, const std::vector< unsigned > &nodes, std::set< std::pair< unsigned, unsigned > > &subCells)
void getInterpolatoryCoordinates(Kokkos::DynRankView< double, PHX::Device > &coords) const
virtual int getSubcellCount(int dim) const
Teuchos::RCP< Intrepid2::Basis< PHX::Device, double, double > > intrepidBasis_