Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GlobalIndexer_Utilities.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 "PanzerDofMgr_config.hpp"
13 
14 namespace panzer {
15 
16 
17 std::vector<Teuchos::RCP<const GlobalIndexer>>
18 nc2c_vector(const std::vector<Teuchos::RCP<GlobalIndexer > > & ugis)
19 {
20  std::vector<Teuchos::RCP<const GlobalIndexer>> vec;
21 
22  for(std::size_t blk=0;blk<ugis.size();blk++)
23  vec.push_back(ugis[blk]);
24 
25  return vec;
26 }
27 
28 int getFieldBlock(const std::string & fieldName,
29  const std::vector<Teuchos::RCP<const GlobalIndexer>> & ugis)
30 {
31  int fieldNum = -1;
32  for(std::size_t blk=0;blk<ugis.size();blk++) {
33  fieldNum = ugis[blk]->getFieldNum(fieldName);
34  if(fieldNum>=0)
35  return blk;
36  }
37 
38  return fieldNum;
39 }
40 
41 int getFieldBlock(const std::string & fieldName,
42  const std::vector<Teuchos::RCP<GlobalIndexer>> & ugis)
43 {
44  int fieldNum = -1;
45  for(std::size_t blk=0;blk<ugis.size();blk++) {
46  fieldNum = ugis[blk]->getFieldNum(fieldName);
47  if(fieldNum>=0)
48  return fieldNum;
49  }
50 
51  return fieldNum;
52 }
53 
54 void computeBlockOffsets(const std::string & blockId,
55  const std::vector<Teuchos::RCP<GlobalIndexer>> & ugis,
56  std::vector<int> & blockOffsets)
57 {
58  blockOffsets.resize(ugis.size()+1); // number of fields, plus a sentinnel
59 
60  int offset = 0;
61  for(std::size_t blk=0;blk<ugis.size();blk++) {
62  blockOffsets[blk] = offset;
63  offset += ugis[blk]->getElementBlockGIDCount(blockId);
64  }
65  blockOffsets[ugis.size()] = offset;
66 }
67 
68 void computeBlockOffsets(const std::string & blockId,
69  const std::vector<Teuchos::RCP<const GlobalIndexer>> & ugis,
70  std::vector<int> & blockOffsets)
71 {
72  blockOffsets.resize(ugis.size()+1); // number of fields, plus a sentinnel
73 
74  int offset = 0;
75  for(std::size_t blk=0;blk<ugis.size();blk++) {
76  blockOffsets[blk] = offset;
77  offset += ugis[blk]->getElementBlockGIDCount(blockId);
78  }
79  blockOffsets[ugis.size()] = offset;
80 }
81 
82 std::string
84 {
85  std::size_t myOwnedCount = static_cast<std::size_t>(ugi.getNumOwned());
86  std::size_t sum=0,min=0,max=0;
87 
88  // get min,max and sum
89  Teuchos::reduceAll(*ugi.getComm(),Teuchos::REDUCE_SUM,1,&myOwnedCount,&sum);
90  Teuchos::reduceAll(*ugi.getComm(),Teuchos::REDUCE_MIN,1,&myOwnedCount,&min);
91  Teuchos::reduceAll(*ugi.getComm(),Teuchos::REDUCE_MAX,1,&myOwnedCount,&max);
92 
93  // compute mean and variance
94  double dev2 = (double(myOwnedCount)-double(sum)/double(ugi.getComm()->getSize()));
95  dev2 *= dev2;
96 
97  double variance = 0.0;
98  Teuchos::reduceAll(*ugi.getComm(),Teuchos::REDUCE_SUM,1,&dev2,&variance);
99 
100  double mean = sum / double(ugi.getComm()->getSize());
101  variance = variance / double(ugi.getComm()->getSize());
102 
103  // now print to a string stream
104  std::stringstream ss;
105  ss << "Max, Min, Mean, StdDev = " << max << ", " << min << ", " << mean << ", " << std::sqrt(variance);
106 
107  return ss.str();
108 }
109 
110 void
111 printMeshTopology(std::ostream & os,const panzer::GlobalIndexer & ugi)
112 {
113  std::vector<std::string> block_ids;
114 
115  ugi.getElementBlockIds(block_ids);
116  for(std::size_t b=0;b<block_ids.size();b++) {
117  // extract the elemnts in each element block
118  const std::vector<panzer::LocalOrdinal> & elements = ugi.getElementBlock(block_ids[b]);
119 
120  os << "Element Block: \"" << block_ids[b] << "\"" << std::endl;
121 
122  // loop over element in this element block, write out to
123  for(std::size_t e=0;e<elements.size();e++) {
124  // extract LIDs, this is terribly inefficient for certain
125  // devices but only used for debugging
126  PHX::View<const int*> lids = ugi.getElementLIDs(elements[e]);
127  auto lids_host = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(),lids);
128 
129  // extract GIDs, this array is filled
130  std::vector<panzer::GlobalOrdinal> gids;
131  ugi.getElementGIDs(elements[e],gids);
132 
133  os << " local element id = " << elements[e] << ", ";
134 
135  os << " gids =";
136  for(std::size_t i=0;i<gids.size();i++)
137  os << " " << gids[i];
138 
139  os << ", lids =";
140  for(std::size_t i=0;i<gids.size();i++)
141  os << " " << lids_host(i);
142  os << std::endl;
143  }
144  }
145 }
146 
149 {
150  using Map = Tpetra::Map<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
151  using IntVector = Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
152 
153  std::vector<panzer::GlobalOrdinal> indices;
154  std::vector<std::string> blocks;
155 
156  ugi.getOwnedAndGhostedIndices(indices);
157  ugi.getElementBlockIds(blocks);
158 
159  std::vector<int> fieldNumbers(indices.size(),-1);
160 
161  Teuchos::RCP<Map> ghostedMap
162  = Teuchos::rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(), Teuchos::arrayViewFromVector(indices),
164 
165  // build a map from local ids to a field number
166  for(std::size_t blk=0;blk<blocks.size();blk++) {
167  std::string blockId = blocks[blk];
168 
169  const std::vector<panzer::LocalOrdinal> & elements = ugi.getElementBlock(blockId);
170  const std::vector<int> & fields = ugi.getBlockFieldNumbers(blockId);
171 
172  // loop over all elements, and set field number in output array
173  std::vector<panzer::GlobalOrdinal> gids(fields.size());
174  for(std::size_t e=0;e<elements.size();e++) {
175  ugi.getElementGIDs(elements[e],gids);
176 
177  for(std::size_t f=0;f<fields.size();f++) {
178  int fieldNum = fields[f];
179  panzer::GlobalOrdinal gid = gids[f];
180  std::size_t lid = ghostedMap->getLocalElement(gid); // hash table lookup
181 
182  fieldNumbers[lid] = fieldNum;
183  }
184  }
185  }
186 
187  // produce a reduced vector containing only fields known by this processor
188  std::vector<panzer::GlobalOrdinal> reducedIndices;
189  std::vector<int> reducedFieldNumbers;
190  for(std::size_t i=0;i<fieldNumbers.size();i++) {
191  if(fieldNumbers[i]>-1) {
192  reducedIndices.push_back(indices[i]);
193  reducedFieldNumbers.push_back(fieldNumbers[i]);
194  }
195  }
196 
197  Teuchos::RCP<Map> reducedMap
198  = Teuchos::rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(), Teuchos::arrayViewFromVector(reducedIndices),
200  return Teuchos::rcp(new IntVector(reducedMap,Teuchos::arrayViewFromVector(reducedFieldNumbers)));
201 }
202 
204  std::vector<int> & fieldNumbers,
205  const Teuchos::RCP<const Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> > & reducedVec)
206 {
207  using IntVector = Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
209 
210  auto host_values = dest->getLocalViewHost(Tpetra::Access::ReadOnly);
211 
212  fieldNumbers.resize(dest->getLocalLength());
213  for (size_t i=0; i < fieldNumbers.size(); ++i)
214  fieldNumbers[i] = host_values(i,0);
215 }
216 
219  const Teuchos::RCP<const Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> > & reducedVec)
220 {
221  using Map = Tpetra::Map<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
222  using IntVector = Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
223  using Importer = Tpetra::Import<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>;
224 
225  // first step: get a reduced field number vector and build a map to
226  // contain the full field number vector
228 
229  Teuchos::RCP<Map> destMap;
230  {
231  std::vector<panzer::GlobalOrdinal> indices;
232  ugi.getOwnedAndGhostedIndices(indices);
233  destMap = Teuchos::rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(), Teuchos::arrayViewFromVector(indices),
235  }
236 
237  Teuchos::RCP<const IntVector> source = reducedVec;
238  if(source==Teuchos::null)
239  source = buildGhostedFieldReducedVector(ugi);
240  Teuchos::RCP<const Map> sourceMap = source->getMap();
241 
242  // second step: perform the global communciation required to fix the
243  // interface conditions (where this processor doesn't know what field
244  // some indices are)
246  Teuchos::RCP<IntVector> dest = Teuchos::rcp(new IntVector(destMap));
247  Importer importer(sourceMap,destMap);
248 
249  dest->doImport(*source,importer,Tpetra::INSERT);
250  PHX::Device::execution_space().fence();
251 
252  return dest;
253 }
254 
257 getFieldMap(int fieldNum,const Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> & fieldTVector)
258 {
260  std::vector<int> fieldVector(fieldTVector.getLocalLength());
261  fieldTVector.get1dCopy(Teuchos::arrayViewFromVector(fieldVector));
262 
263  std::vector<panzer::GlobalOrdinal> mapVector;
264  for(std::size_t i=0;i<fieldVector.size();i++) {
265  if(fieldVector[i]==fieldNum)
266  mapVector.push_back(origMap->getGlobalElement(i));
267  }
268 
270  = Teuchos::rcp(new Tpetra::Map<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(
271  Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(), Teuchos::arrayViewFromVector(mapVector),
273 
274  return finalMap;
275 }
276 
278  : ugi_(ugi)
279 {
282 }
283 
284 
285 void ArrayToFieldVector::buildFieldVector(const Tpetra::Vector<int,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> & source) const
286 {
287  // build (unghosted) vector and map
288  std::vector<panzer::GlobalOrdinal> indices;
289  ugi_->getOwnedIndices(indices);
290 
292  = Teuchos::rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(), Teuchos::arrayViewFromVector(indices),
294  Teuchos::RCP<IntVector> localFieldVector = Teuchos::rcp(new IntVector(destMap));
295 
296  Tpetra::Import<int,panzer::GlobalOrdinal,panzer::TpetraNodeType> importer(source.getMap(),destMap);
297  localFieldVector->doImport(source,importer,Tpetra::INSERT);
298 
299  fieldVector_ = localFieldVector;
300 }
301 
303 ArrayToFieldVector::getFieldMap(const std::string & fieldName) const
304 {
305  return getFieldMap(ugi_->getFieldNum(fieldName));
306 }
307 
310 {
311  if(fieldMaps_[fieldNum]==Teuchos::null) {
312  // if neccessary build field vector
313  if(fieldVector_==Teuchos::null)
315 
316  fieldMaps_[fieldNum] = panzer::getFieldMap(fieldNum,*fieldVector_);
317  }
318 
319  return fieldMaps_[fieldNum];
320 }
321 
322 // *****************************************
323 // namespace orientation_helpers
324 // *****************************************
325 
326 namespace orientation_helpers {
327 
328 
329 void computeCellEdgeOrientations(const std::vector<std::pair<int,int> > & topEdgeIndices,
330  const std::vector<panzer::GlobalOrdinal> & topology,
331  const FieldPattern & fieldPattern,
332  std::vector<signed char> & orientation)
333 {
334  // LOCAL element orientations are always set so that they flow in the positive
335  // direction along an edge from node 0 to node 1. As a result if the GID of
336  // node 0 is larger then node 1 then the GLOBAL orientation is -1 (and positive
337  // otherwise). The local definition of the edge direction is defined by
338  // the shards cell topology.
339 
340  TEUCHOS_ASSERT(orientation.size()==std::size_t(fieldPattern.numberIds()));
341 
342  int edgeDim = 1;
343 
344  for(std::size_t e=0;e<topEdgeIndices.size();e++) {
345  // grab topological nodes
346  const std::pair<int,int> nodes = topEdgeIndices[e];
347 
348  // extract global values of topological nodes
349  panzer::GlobalOrdinal v0 = topology[nodes.first];
350  panzer::GlobalOrdinal v1 = topology[nodes.second];
351 
352  // using simple rule make a decision about orientation
353  signed char edgeOrientation = 1;
354  if(v1>v0)
355  edgeOrientation = 1;
356  else if(v0>v1)
357  edgeOrientation = -1;
358  else
359  { TEUCHOS_ASSERT(false); }
360 
361  // grab edgeIndices to be set to compute orientation
362  const std::vector<int> & edgeIndices = fieldPattern.getSubcellIndices(edgeDim,e);
363  for(std::size_t s=0;s<edgeIndices.size();s++)
364  orientation[edgeIndices[s]] = edgeOrientation;
365  }
366 }
367 
368 void computeCellFaceOrientations(const std::vector<std::vector<int> > & topFaceIndices,
369  const std::vector<panzer::GlobalOrdinal> & topology,
370  const FieldPattern & fieldPattern,
371  std::vector<signed char> & orientation)
372 {
373  // LOCAL element orientations are always set so that they flow in the positive
374  // direction away from the cell (counter clockwise rotation of the face). To determine
375  // the face orientation we use the fact that Shards (and thus the field pattern) always
376  // locally orders faces in a counter clockwise direction. A local rule for each element
377  // will take advantage of this ensuring both elements agree on the orientation. This rule
378  // is to first find the smallest node GID on the face. Then look at the GIDs of the nodes
379  // immediately preceding and following that one. If the node following has smaller GID than
380  // the preceding node then the face is oriented counter clockwise and thus the orientation
381  // is +1. If the node preceding is larger, then the orientation is clockwise and the set to
382  // a value of -1.
383 
384  // this only works for 3D field patterns
385  TEUCHOS_ASSERT(fieldPattern.getDimension()==3);
386 
387  TEUCHOS_ASSERT(orientation.size()==std::size_t(fieldPattern.numberIds()));
388 
389  int faceDim = 2;
390 
391  for(std::size_t f=0;f<topFaceIndices.size();f++) {
392  // grab topological nodes
393  const std::vector<int> & nodes = topFaceIndices[f];
394  std::vector<panzer::GlobalOrdinal> globals(nodes.size());
395  for(std::size_t n=0;n<nodes.size();n++)
396  globals[n] = topology[nodes[n]];
397 
398  typename std::vector<panzer::GlobalOrdinal>::const_iterator itr
399  = std::min_element(globals.begin(),globals.end());
400 
401  TEUCHOS_TEST_FOR_EXCEPTION(itr==globals.end(),std::out_of_range,
402  "panzer::orientation_helpers::computeCellFaceOrientations: A face index array "
403  "was empty.");
404 
405  // extract global values of topological nodes
406  // The face nodes go in counter clockwise order, let v_min be the
407  // value with the minimum element then
408  // vbefore => itr => vafter
409  // note that the nonsense with the beginning and end has to do with
410  // if this iterator was the first or last in the array
411  panzer::GlobalOrdinal vbefore = itr==globals.begin() ? *(globals.end()-1) : *(itr-1);
412  panzer::GlobalOrdinal vafter = (itr+1)==globals.end() ? *globals.begin() : *(itr+1);
413 
414 /*
415  // sanity check in debug mode (uncomment these lines)
416  TEUCHOS_ASSERT(std::find(globals.begin(),globals.end(),vbefore)!=globals.end());
417  TEUCHOS_ASSERT(std::find(globals.begin(),globals.end(),vafter)!=globals.end());
418 
419  // print out information about the found nodes and also what
420  // order they were in originally
421  std::cout << "\nFace Order = ";
422  for(std::size_t l=0;l<globals.size();l++)
423  std::cout << globals[l] << " ";
424  std::cout << std::endl;
425  std::cout << "(before,min,after) " << f << ": " << vbefore << " => " << *itr << " => " << vafter << std::endl;
426 */
427 
428  // note by assumption
429  // vbefore < *itr and *itr < vafter
430 
431  // Based on the next lowest global id starting from the minimum
432  signed char faceOrientation = 1;
433  if(vafter>vbefore) // means smaller in clockwise direction
434  faceOrientation = -1;
435  else if(vbefore>vafter) // means smaller in counter clockwise direction
436  faceOrientation = 1;
437  else
438  { TEUCHOS_ASSERT(false); } // we got an equality somehow!
439 
440  // grab faceIndices to be set to compute orientation
441  const std::vector<int> & faceIndices = fieldPattern.getSubcellIndices(faceDim,f);
442  for(std::size_t s=0;s<faceIndices.size();s++)
443  orientation[faceIndices[s]] = faceOrientation;
444  }
445 }
446 
447 void computePatternEdgeIndices(const FieldPattern & pattern,std::vector<std::pair<int,int> > & edgeIndices)
448 {
449  unsigned dim = 1;
450  shards::CellTopology cellTopo = pattern.getCellTopology();
451  for(unsigned e=0;e<cellTopo.getEdgeCount();e++) {
452  // get local vertex ids for a this edge
453  unsigned local_v0 = cellTopo.getNodeMap(dim,e,0);
454  unsigned local_v1 = cellTopo.getNodeMap(dim,e,1);
455 
456  // get sub cell indices for geometric pattern
457  const std::vector<int> & v0_indices = pattern.getSubcellIndices(0,local_v0);
458  const std::vector<int> & v1_indices = pattern.getSubcellIndices(0,local_v1);
459 
460  TEUCHOS_ASSERT(v0_indices.size()>0); // there must be a node
461  TEUCHOS_ASSERT(v1_indices.size()>0); // there must be a node
462 
463  // take the first index on each vertex and make a edge lookup
464  edgeIndices.push_back(std::make_pair(v0_indices[0],v1_indices[0]));
465  }
466 }
467 
468 void computePatternFaceIndices(const FieldPattern & pattern,std::vector<std::vector<int> > & faceIndices)
469 {
470  // this only works for 3D field patterns
471  // TEUCHOS_ASSERT(pattern.getDimension()==3);
472  //
473  unsigned node_dim = 0; // by assumption
474  unsigned subcell_dim = 2;
475 
476  if(pattern.getDimension()==3) {
477  shards::CellTopology cellTopo = pattern.getCellTopology();
478 
479  faceIndices.resize(cellTopo.getSubcellCount(subcell_dim));
480 
481  for(unsigned f=0;f<cellTopo.getSubcellCount(subcell_dim);f++) {
482  shards::CellTopology faceTopo(cellTopo.getBaseCellTopologyData(subcell_dim,f));
483 
484  for(unsigned v=0;v<faceTopo.getNodeCount();v++) {
485  // get local vertex ids for a this edge
486  unsigned local_v = cellTopo.getNodeMap(subcell_dim,f,v);
487 
488  // get sub cell indices for geometric pattern
489  const std::vector<int> & v_indices = pattern.getSubcellIndices(node_dim,local_v);
490 
491  TEUCHOS_ASSERT(v_indices.size()>0); // there must be a node
492 
493  // take the first index on each vertex and make a edge lookup
494  faceIndices[f].push_back(v_indices[0]);
495  }
496  }
497  }
498  else if(pattern.getDimension()==2) {
499  shards::CellTopology cellTopo = pattern.getCellTopology();
500 
501  faceIndices.resize(1);
502 
503  for(unsigned v=0;v<cellTopo.getNodeCount();v++)
504  faceIndices[0].push_back(v);
505  }
506 }
507 
508 } // end namespace orientation_helpers
509 } // end namespace panzer
std::map< int, Teuchos::RCP< const Map > > fieldMaps_
(unghosted) field vector (as needed)
virtual Teuchos::RCP< Teuchos::Comm< int > > getComm() const =0
virtual int getDimension() const =0
const Kokkos::View< const panzer::LocalOrdinal *, Kokkos::LayoutRight, PHX::Device > getElementLIDs(panzer::LocalOrdinal localElmtId) const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const IntVector > gh_reducedFieldVector_
virtual const std::vector< panzer::LocalOrdinal > & getElementBlock(const std::string &blockId) const =0
Teuchos::RCP< Tpetra::Vector< int, int, panzer::GlobalOrdinal, panzer::TpetraNodeType > > buildGhostedFieldReducedVector(const GlobalIndexer &ugi)
virtual void getElementBlockIds(std::vector< std::string > &elementBlockIds) const =0
Teuchos::RCP< const IntVector > gh_fieldVector_
ghosted reduced field vector
Teuchos::RCP< const Tpetra::Map< int, panzer::GlobalOrdinal, panzer::TpetraNodeType > > getFieldMap(int fieldNum, const Tpetra::Vector< int, int, panzer::GlobalOrdinal, panzer::TpetraNodeType > &fieldTVector)
void computeBlockOffsets(const std::string &blockId, const std::vector< Teuchos::RCP< GlobalIndexer >> &ugis, std::vector< int > &blockOffsets)
virtual int numberIds() const
Tpetra::Vector< int, int, panzer::GlobalOrdinal, panzer::TpetraNodeType > IntVector
Kokkos::View< const LO **, Kokkos::LayoutRight, PHX::Device > lids
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< const GlobalIndexer > ugi_
DOF mapping.
std::vector< Teuchos::RCP< const GlobalIndexer > > nc2c_vector(const std::vector< Teuchos::RCP< GlobalIndexer > > &ugis)
int getFieldBlock(const std::string &fieldName, const std::vector< Teuchos::RCP< const GlobalIndexer >> &ugis)
void computeCellEdgeOrientations(const std::vector< std::pair< int, int > > &topEdgeIndices, const std::vector< panzer::GlobalOrdinal > &topology, const FieldPattern &fieldPattern, std::vector< signed char > &orientation)
virtual void getOwnedAndGhostedIndices(std::vector< panzer::GlobalOrdinal > &indices) const =0
Get the set of owned and ghosted indices for this processor.
void computePatternFaceIndices(const FieldPattern &pattern, std::vector< std::vector< int > > &faceIndices)
virtual int getNumOwned() const =0
Get the number of indices owned by this processor.
void computeCellFaceOrientations(const std::vector< std::vector< int > > &topFaceIndices, const std::vector< panzer::GlobalOrdinal > &topology, const FieldPattern &fieldPattern, std::vector< signed char > &orientation)
Tpetra::Map< int, panzer::GlobalOrdinal, panzer::TpetraNodeType > Map
virtual shards::CellTopology getCellTopology() const =0
Teuchos::RCP< const Tpetra::Map< int, panzer::GlobalOrdinal, panzer::TpetraNodeType > > getFieldMap(const std::string &fieldName) const
virtual const std::vector< int > & getBlockFieldNumbers(const std::string &blockId) const =0
virtual void getElementGIDs(panzer::LocalOrdinal localElmtId, std::vector< panzer::GlobalOrdinal > &gids, const std::string &blockIdHint="") const =0
Get the global IDs for a particular element. This function overwrites the gids variable.
virtual const std::vector< int > & getSubcellIndices(int dim, int cellIndex) const =0
void buildFieldVector(const Tpetra::Vector< int, int, panzer::GlobalOrdinal, panzer::TpetraNodeType > &source) const
build unghosted field vector from ghosted field vector
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< const IntVector > fieldVector_
Maps for each field (as needed)
void printMeshTopology(std::ostream &os, const panzer::GlobalIndexer &ugi)
void computePatternEdgeIndices(const FieldPattern &pattern, std::vector< std::pair< int, int > > &edgeIndices)
void buildGhostedFieldVector(const GlobalIndexer &ugi, std::vector< int > &fieldNumbers, const Teuchos::RCP< const Tpetra::Vector< int, int, panzer::GlobalOrdinal, panzer::TpetraNodeType > > &reducedVec)
ArrayToFieldVector()
Maps for each field (as needed)
std::string printUGILoadBalancingInformation(const GlobalIndexer &ugi)