Panzer Version of the Day
Loading...
Searching...
No Matches
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
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
25namespace panzer {
26
31 public:
32 Intrepid2FieldPattern(const Teuchos::RCP< Intrepid2::Basis<PHX::Device,double,double> > &intrepidBasis);
33 Teuchos::RCP<panzer::FieldPattern> clone() const
34 {return Teuchos::rcp(new Intrepid2FieldPattern(*this));}
35
36 virtual int getSubcellCount(int dim) const;
37 virtual const std::vector<int> & getSubcellIndices(int dim, int cellIndex) const;
38 virtual int getDimension() const;
39 virtual shards::CellTopology getCellTopology() const;
40
41 virtual void getSubcellClosureIndices(int dim,int cellIndex,std::vector<int> & indices) const;
42
43 // static functions for examining shards objects
44
56 static void buildSubcellClosure(const shards::CellTopology & cellTopo,unsigned dim,unsigned subCell,
57 std::set<std::pair<unsigned,unsigned> > & closure);
58
70 static void findContainedSubcells(const shards::CellTopology & cellTopo,unsigned dim,
71 const std::vector<unsigned> & nodes,
72 std::set<std::pair<unsigned,unsigned> > & subCells);
73
81 static void getSubcellNodes(const shards::CellTopology & cellTopo,unsigned dim,unsigned subCell,
82 std::vector<unsigned> & nodes);
83
92
98 void getInterpolatoryCoordinates(Kokkos::DynRankView<double,PHX::Device> & coords) const;
99
109 void getInterpolatoryCoordinates(const Kokkos::DynRankView<double,PHX::Device> & cellNodes,
110 Kokkos::DynRankView<double,PHX::Device> & coords,
111 Teuchos::RCP<const shards::CellTopology> meshCellTopology=Teuchos::null) const;
112
114 Teuchos::RCP< Intrepid2::Basis<PHX::Device,double,double> > getIntrepidBasis() const;
115
116 protected:
117 Teuchos::RCP< Intrepid2::Basis<PHX::Device,double,double> > intrepidBasis_;
118
119 //mutable std::vector<int> subcellIndices_;
120 mutable std::vector<std::vector<std::vector<int> > > subcellIndicies_;
121 std::vector<int> empty_;
122 };
123
124}
125
126#endif
Teuchos::RCP< Intrepid2::Basis< PHX::Device, double, double > > intrepidBasis_
virtual void getSubcellClosureIndices(int dim, int cellIndex, std::vector< int > &indices) const
static void getSubcellNodes(const shards::CellTopology &cellTopo, unsigned dim, unsigned subCell, std::vector< unsigned > &nodes)
void getInterpolatoryCoordinates(Kokkos::DynRankView< double, PHX::Device > &coords) const
virtual shards::CellTopology getCellTopology() const
Teuchos::RCP< Intrepid2::Basis< PHX::Device, double, double > > getIntrepidBasis() const
Returns the underlying Intrepid2::Basis object.
bool supportsInterpolatoryCoordinates() const
Does this field pattern support interpolatory coordinates?
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)
virtual const std::vector< int > & getSubcellIndices(int dim, int cellIndex) const
Teuchos::RCP< panzer::FieldPattern > clone() const
static void buildSubcellClosure(const shards::CellTopology &cellTopo, unsigned dim, unsigned subCell, std::set< std::pair< unsigned, unsigned > > &closure)