Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_IntrepidOrientation.cpp
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 #include "PanzerDiscFE_config.hpp"
12 #include "Panzer_GlobalIndexer.hpp"
14 
15 namespace panzer {
16 
17  void buildIntrepidOrientation(std::vector<Intrepid2::Orientation> & orientation,
18  panzer::ConnManager & connMgr)
19  {
20  using Teuchos::rcp_dynamic_cast;
21  using Teuchos::RCP;
22  using Teuchos::rcp;
23 
24  orientation.clear();
25 
26  // Retrive element blocks and its meta data
27  const int numElementBlocks = connMgr.numElementBlocks();
28 
29  std::vector<std::string> elementBlockIds;
30  std::vector<shards::CellTopology> elementBlockTopologies;
31 
32  connMgr.getElementBlockIds(elementBlockIds);
33  connMgr.getElementBlockTopologies(elementBlockTopologies);
34 
35  TEUCHOS_TEST_FOR_EXCEPTION(numElementBlocks <= 0 &&
36  numElementBlocks != static_cast<int>(elementBlockIds.size()) &&
37  numElementBlocks != static_cast<int>(elementBlockTopologies.size()),
38  std::logic_error,
39  "panzer::buildIntrepidOrientation: Number of element blocks does not match to element block meta data");
40 
41  // Currently panzer support only one type of elements for whole mesh (use the first cell topology)
42  const auto cellTopo = elementBlockTopologies.at(0);
43  const int numVerticesPerCell = cellTopo.getVertexCount();
44 
45  const auto fp = NodalFieldPattern(cellTopo);
46  connMgr.buildConnectivity(fp);
47 
48  // Count and pre-alloc orientations
49  int total_elems = 0;
50  for (int i=0;i<numElementBlocks;++i) {
51  total_elems += connMgr.getElementBlock(elementBlockIds.at(i)).size();
52  }
53 
54  orientation.resize(total_elems);
55  // Loop over element blocks
56  for (int i=0;i<numElementBlocks;++i) {
57  // get elements in a block
58  const auto &elementBlock = connMgr.getElementBlock(elementBlockIds.at(i));
59 
60  const int numElementsPerBlock = elementBlock.size();
61 
62  // construct orientation information
63  for (int c=0;c<numElementsPerBlock;++c) {
64  const int localCellId = elementBlock.at(c);
66  vertices(connMgr.getConnectivity(localCellId), numVerticesPerCell);
67  // This function call expects a view for the vertices, not the nodes
68  orientation[localCellId] = (Intrepid2::Orientation::getOrientation(cellTopo, vertices));
69  }
70  }
71  }
72 
75  {
76  using Teuchos::rcp_dynamic_cast;
77  using Teuchos::RCP;
78  using Teuchos::rcp;
79 
80  auto orientation = rcp(new std::vector<Intrepid2::Orientation>);
81 
82  {
84  = rcp_dynamic_cast<const GlobalIndexer>(globalIndexer);
85 
86  if (ugi!=Teuchos::null) {
87  const auto connMgr = ugi->getConnManager()->noConnectivityClone();
88 
89  TEUCHOS_TEST_FOR_EXCEPTION(connMgr == Teuchos::null,std::logic_error,
90  "panzer::buildIntrepidOrientation: ConnManager is null!");
91 
92  buildIntrepidOrientation(*orientation, *connMgr);
93  return orientation;
94  }
95  }
96 
97  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
98  "panzer::buildIntrepidOrientation: Could not cast GlobalIndexer");
99  }
100 
101  void buildIntrepidOrientations(const std::vector<std::string>& eBlockNames,
102  const panzer::ConnManager & connMgrInput,
103  std::map<std::string,std::vector<Intrepid2::Orientation>> & orientations)
104  {
105  using Teuchos::rcp_dynamic_cast;
106  using Teuchos::RCP;
107  using Teuchos::rcp;
108 
109  auto connMgrPtr = connMgrInput.noConnectivityClone();
110  auto& connMgr = *connMgrPtr;
111 
112  // Map element block name to topology
113  std::map<std::string,shards::CellTopology> eb_name_to_topo;
114  {
115  std::vector<std::string> elementBlockIds;
116  connMgr.getElementBlockIds(elementBlockIds);
117 
118  std::vector<shards::CellTopology> elementBlockTopologies;
119  connMgr.getElementBlockTopologies(elementBlockTopologies);
120 
121  for (size_t i=0; i < elementBlockIds.size(); ++i) {
122  eb_name_to_topo[elementBlockIds[i]] = elementBlockTopologies[i];
123  }
124  }
125 
126  // Currently panzer supports only one element topology for whole mesh (use the first cell topology)
127  const auto cellTopo = eb_name_to_topo[eBlockNames[0]];
128  const int numVerticesPerCell = cellTopo.getVertexCount();
129 
130  const auto fp = NodalFieldPattern(cellTopo);
131  connMgr.buildConnectivity(fp);
132 
133  // Loop over requested element blocks
134  for (size_t i=0;i<eBlockNames.size();++i) {
135 
136  const auto &lids = connMgr.getElementBlock(eBlockNames[i]);
137  const size_t num_elems = lids.size();
138  auto& orts = orientations[eBlockNames[i]];
139  orts.resize(num_elems);
140 
141  for (size_t c=0;c<num_elems;++c) {
142  const int lid = lids[c];
144  vertices(connMgr.getConnectivity(lid),numVerticesPerCell);
145 
146  // This function call expects a view for the vertices, not the nodes
147  orts[c] = Intrepid2::Orientation::getOrientation(cellTopo, vertices);
148  }
149  }
150  }
151 
152 } // end namespace panzer
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void buildIntrepidOrientations(const std::vector< std::string > &eBlockNames, const panzer::ConnManager &connMgrInput, std::map< std::string, std::vector< Intrepid2::Orientation >> &orientations)
Builds the element orientations for the specified element blocks.
void buildIntrepidOrientation(std::vector< Intrepid2::Orientation > &orientation, panzer::ConnManager &connMgr)
Builds the element orientations for all element blocks.
Kokkos::View< const LO **, Kokkos::LayoutRight, PHX::Device > lids
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual const std::vector< LocalOrdinal > & getElementBlock(const std::string &blockID) const =0
virtual void getElementBlockIds(std::vector< std::string > &elementBlockIds) const =0
virtual Teuchos::RCP< ConnManager > noConnectivityClone() const =0
virtual Teuchos::RCP< const ConnManager > getConnManager() const =0
Returns the connection manager currently being used.
Pure virtual base class for supplying mesh connectivity information to the DOF Manager.
virtual void getElementBlockTopologies(std::vector< shards::CellTopology > &elementBlockTopologies) const =0
virtual std::size_t numElementBlocks() const =0
virtual void buildConnectivity(const FieldPattern &fp)=0
virtual const GlobalOrdinal * getConnectivity(LocalOrdinal localElmtId) const =0