Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_SetupUtilities.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
12#include "Panzer_Workset_Builder.hpp"
13#include "Teuchos_Assert.hpp"
14
15#include <stk_mesh/base/Selector.hpp>
16#include <stk_mesh/base/GetEntities.hpp>
17
18namespace panzer_stk {
19
20Teuchos::RCP<std::vector<panzer::Workset> >
22 const std::string & eBlock,
23 const panzer::WorksetNeeds & needs)
24{
25 using namespace workset_utils;
26
27 std::vector<std::string> element_blocks;
28
29 std::vector<std::size_t> local_cell_ids;
30 Kokkos::DynRankView<double,PHX::Device> cell_node_coordinates;
31
32 getIdsAndNodes(mesh, eBlock, local_cell_ids, cell_node_coordinates);
33
34 // only build workset if there are elements to worry about
35 // this may be processor dependent, so an element block
36 // may not have elements and thus no contribution
37 // on this processor
38 return panzer::buildWorksets(needs, eBlock, local_cell_ids, cell_node_coordinates);
39}
40
41Teuchos::RCP<std::vector<panzer::Workset> >
43 const panzer::WorksetNeeds & needs,
44 const std::string & sideset,
45 const std::string & eBlock,
46 bool useCascade)
47{
48 using namespace workset_utils;
49 using Teuchos::RCP;
50
51 std::vector<stk::mesh::Entity> sideEntities;
52
53 try {
54 // grab local entities on this side
55 // ...catch any failure...primarily wrong side set and element block info
56 if(!useCascade)
57 mesh.getMySides(sideset,eBlock,sideEntities);
58 else
59 mesh.getAllSides(sideset,eBlock,sideEntities);
60 }
62 std::stringstream ss;
63 std::vector<std::string> sideSets;
64 mesh.getSidesetNames(sideSets);
65
66 // build an error message
67 ss << e.what() << "\nChoose one of:\n";
68 for(std::size_t i=0;i<sideSets.size();i++)
69 ss << "\"" << sideSets[i] << "\"\n";
70
71 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
72 }
74 std::stringstream ss;
75 std::vector<std::string> elementBlocks;
76 mesh.getElementBlockNames(elementBlocks);
77
78 // build an error message
79 ss << e.what() << "\nChoose one of:\n";
80 for(std::size_t i=0;i<elementBlocks.size();i++)
81 ss << "\"" << elementBlocks[i] << "\"\n";
82
83 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
84 }
85 catch(std::logic_error & e) {
86 std::stringstream ss;
87 ss << e.what() << "\nUnrecognized logic error.\n";
88
89 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
90 }
91
92 std::vector<stk::mesh::Entity> elements;
93 std::map<std::pair<unsigned,unsigned>,std::vector<std::size_t> > local_cell_ids;
94 if(!useCascade) {
95 unsigned subcell_dim = needs.cellData.baseCellDimension()-1;
96 std::vector<std::size_t> local_side_ids;
97 getSideElements(mesh, eBlock,
98 sideEntities,local_side_ids,elements);
99
100 // build local cell_ids, mapped by local side id
101 for(std::size_t elm=0;elm<elements.size();++elm) {
102 stk::mesh::Entity element = elements[elm];
103
104 local_cell_ids[std::make_pair(subcell_dim,local_side_ids[elm])].push_back(mesh.elementLocalId(element));
105 }
106 }
107 else {
108 std::vector<std::size_t> local_subcell_ids, subcell_dim;
109 getSideElementCascade(mesh, eBlock,
110 sideEntities,subcell_dim,local_subcell_ids,elements);
111
112 // build local cell_ids, mapped by local side id
113 for(std::size_t elm=0;elm<elements.size();++elm) {
114 stk::mesh::Entity element = elements[elm];
115
116 local_cell_ids[std::make_pair(subcell_dim[elm],local_subcell_ids[elm])].push_back(mesh.elementLocalId(element));
117 }
118 }
119
120 // only build workset if there are elements to worry about
121 // this may be processor dependent, so a defined boundary
122 // condition may have not elements and thus no contribution
123 // on this processor
124 if(elements.size()!=0) {
125 Teuchos::RCP<const shards::CellTopology> topo = mesh.getCellTopology(eBlock);
126
127 // worksets to be returned
128 Teuchos::RCP<std::vector<panzer::Workset> > worksets = Teuchos::rcp(new std::vector<panzer::Workset>);
129
130 // loop over each side
131 for(std::map<std::pair<unsigned,unsigned>,std::vector<std::size_t> >::const_iterator itr=local_cell_ids.begin();
132 itr!=local_cell_ids.end();++itr) {
133
134 if(itr->second.size()==0)
135 continue;
136
137 Kokkos::DynRankView<double,PHX::Device> nodes;
138 mesh.getElementNodes(itr->second,eBlock,nodes);
139
140 Teuchos::RCP<std::vector<panzer::Workset> > current
141 = panzer::buildWorksets(needs, eBlock, itr->second, nodes);
142
143 // correct worksets so the sides are correct
144 for(std::size_t w=0;w<current->size();w++) {
145 (*current)[w].subcell_dim = itr->first.first;
146 (*current)[w].subcell_index = itr->first.second;
147 }
148
149 // append new worksets
150 worksets->insert(worksets->end(),current->begin(),current->end());
151 }
152
153 return worksets;
154 }
155
156 // return Teuchos::null;
157 return Teuchos::rcp(new std::vector<panzer::Workset>());
158}
159
160Teuchos::RCP<std::map<unsigned,panzer::Workset> >
162 const panzer::WorksetNeeds & needs_a,
163 const std::string & blockid_a,
164 const panzer::WorksetNeeds & needs_b,
165 const std::string & blockid_b,
166 const std::string & sideset)
167{
168 using namespace workset_utils;
169 using Teuchos::RCP;
170
171 std::vector<stk::mesh::Entity> sideEntities; // we will reduce a_ and b_ to this vector
172
173 try {
174 // grab local entities on this side
175 // ...catch any failure...primarily wrong side set and element block info
176
177 // we can't use getMySides because it only returns locally owned sides
178 // this gurantees all the sides are extracted (element ownership is considered
179 // we we call getSideElements below)
180
181 stk::mesh::Part * sidePart = mesh.getSideset(sideset);
182 TEUCHOS_TEST_FOR_EXCEPTION(sidePart==0,std::logic_error,
183 "Unknown side set \"" << sideset << "\"");
184
185 stk::mesh::Selector side = *sidePart;
186 // stk::mesh::Selector ownedBlock = metaData_->locally_owned_part() & side;
187
188 // grab elements
189 stk::mesh::get_selected_entities(side,mesh.getBulkData()->buckets(mesh.getSideRank()),sideEntities);
190 }
192 std::stringstream ss;
193 std::vector<std::string> elementBlocks;
194 mesh.getElementBlockNames(elementBlocks);
195
196 // build an error message
197 ss << e.what() << "\nChoose one of:\n";
198 for(std::size_t i=0;i<elementBlocks.size();i++)
199 ss << "\"" << elementBlocks[i] << "\"\n";
200
201 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
202 }
203 catch(std::logic_error & e) {
204 std::stringstream ss;
205 ss << e.what() << "\nUnrecognized logic error.\n";
206
207 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
208 }
209
210 std::vector<stk::mesh::Entity> elements_a, elements_b;
211 std::vector<std::size_t> local_cell_ids_a, local_cell_ids_b;
212 std::vector<std::size_t> local_side_ids_a, local_side_ids_b;
213
214 // this enforces that "a" elements must be owned.
215 getSideElements(mesh, blockid_a,blockid_b, sideEntities,
216 local_side_ids_a,elements_a,
217 local_side_ids_b,elements_b);
218
219 TEUCHOS_TEST_FOR_EXCEPTION(elements_a.size()!=elements_b.size(),std::logic_error,
220 "For a DG type boundary, the number of elements on the \"left\" and \"right\" is not the same.");
221
222 // only build workset if there are elements to worry about
223 // this may be processor dependent, so a defined boundary
224 // condition may have not elements and thus no contribution
225 // on this processor
226 if(elements_a.size()==0)
227 return Teuchos::rcp(new std::map<unsigned,panzer::Workset>);
228
229 // loop over elements of this block (note the assures that element_a and element_b
230 // are the same size, the ordering is the same because the order of sideEntities is
231 // the same
232 for(std::size_t elm=0;elm<elements_a.size();++elm) {
233 stk::mesh::Entity element_a = elements_a[elm];
234 stk::mesh::Entity element_b = elements_b[elm];
235
236 local_cell_ids_a.push_back(mesh.elementLocalId(element_a));
237 local_cell_ids_b.push_back(mesh.elementLocalId(element_b));
238 }
239
240 Kokkos::DynRankView<double,PHX::Device> node_coordinates_a, node_coordinates_b;
241 mesh.getElementNodes(local_cell_ids_a,blockid_a,node_coordinates_a);
242 mesh.getElementNodes(local_cell_ids_b,blockid_b,node_coordinates_b);
243
244 // worksets to be returned
245 return buildBCWorkset(needs_a,blockid_a, local_cell_ids_a, local_side_ids_a, node_coordinates_a,
246 needs_b,blockid_b, local_cell_ids_b, local_side_ids_b, node_coordinates_b);
247}
248
249Teuchos::RCP<std::map<unsigned,panzer::Workset> >
251 const panzer::WorksetNeeds & needs,
252 const std::string & eblockID,
253 const std::string & sidesetID)
254{
255 using namespace workset_utils;
256 using Teuchos::RCP;
257
258 std::vector<stk::mesh::Entity> sideEntities;
259
260 try {
261 // grab local entities on this side
262 // ...catch any failure...primarily wrong side set and element block info
263 mesh.getAllSides(sidesetID,eblockID,sideEntities);
264 }
266 std::stringstream ss;
267 std::vector<std::string> sideSets;
268 mesh.getSidesetNames(sideSets);
269
270 // build an error message
271 ss << e.what() << "\nChoose one of:\n";
272 for(std::size_t i=0;i<sideSets.size();i++)
273 ss << "\"" << sideSets[i] << "\"\n";
274
275 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
276 }
278 std::stringstream ss;
279 std::vector<std::string> elementBlocks;
280 mesh.getElementBlockNames(elementBlocks);
281
282 // build an error message
283 ss << e.what() << "\nChoose one of:\n";
284 for(std::size_t i=0;i<elementBlocks.size();i++)
285 ss << "\"" << elementBlocks[i] << "\"\n";
286
287 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
288 }
289 catch(std::logic_error & e) {
290 std::stringstream ss;
291 ss << e.what() << "\nUnrecognized logic error.\n";
292
293 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
294 }
295
296 std::vector<stk::mesh::Entity> elements;
297 std::vector<std::size_t> local_cell_ids;
298 std::vector<std::size_t> local_side_ids;
299 getSideElements(mesh, eblockID,sideEntities,local_side_ids,elements);
300
301 // loop over elements of this block
302 for(std::size_t elm=0;elm<elements.size();++elm) {
303 stk::mesh::Entity element = elements[elm];
304
305 local_cell_ids.push_back(mesh.elementLocalId(element));
306 }
307
308 // only build workset if there are elements to worry about
309 // this may be processor dependent, so a defined boundary
310 // condition may have not elements and thus no contribution
311 // on this processor
312 if(elements.size()!=0) {
313 Teuchos::RCP<const shards::CellTopology> topo
314 = mesh.getCellTopology(eblockID);
315
316 Kokkos::DynRankView<double,PHX::Device> nodes;
317 mesh.getElementNodes(local_cell_ids,eblockID,nodes);
318
319 return panzer::buildBCWorkset(needs, eblockID, local_cell_ids, local_side_ids, nodes);
320 }
321
322 return Teuchos::null;
323}
324
325namespace workset_utils {
326
328 const std::string & blockId,
329 const std::vector<stk::mesh::Entity> & entities,
330 std::vector<std::size_t> & localEntityIds,
331 std::vector<stk::mesh::Entity> & elements)
332{
333 // for verifying that an element is in specified block
334 stk::mesh::Part * blockPart = mesh.getElementBlockPart(blockId);
335 stk::mesh::Part * ownedPart = mesh.getOwnedPart();
336 stk::mesh::BulkData& bulkData = *mesh.getBulkData();
337
338 // loop over each entitiy extracting elements and local entity ID that
339 // are containted in specified block.
340 std::vector<stk::mesh::Entity>::const_iterator entityItr;
341 for(entityItr=entities.begin();entityItr!=entities.end();++entityItr) {
342 stk::mesh::Entity entity = *entityItr;
343
344 const size_t num_rels = bulkData.num_elements(entity);
345 stk::mesh::Entity const* relations = bulkData.begin_elements(entity);
346 stk::mesh::ConnectivityOrdinal const* ordinals = bulkData.begin_element_ordinals(entity);
347 for(std::size_t e=0; e<num_rels; ++e) {
348 stk::mesh::Entity element = relations[e];
349 std::size_t entityId = ordinals[e];
350
351 // is this element in requested block
352 stk::mesh::Bucket const& bucket = bulkData.bucket(element);
353 bool inBlock = bucket.member(*blockPart);
354 bool onProc = bucket.member(*ownedPart);
355 if(inBlock && onProc) {
356 // add element and Side ID to output vectors
357 elements.push_back(element);
358 localEntityIds.push_back(entityId);
359 }
360 }
361 }
362}
363
365 const std::string & blockId,
366 const std::vector<stk::mesh::Entity> & entities,
367 std::vector<std::size_t> & localEntityIds,
368 std::vector<stk::mesh::Entity> & elements,
369 std::vector<std::size_t> & missingElementIndices)
370{
371 // for verifying that an element is in specified block
372 stk::mesh::Part * blockPart = mesh.getElementBlockPart(blockId);
373 stk::mesh::Part * universalPart = &mesh.getMetaData()->universal_part();
374 stk::mesh::BulkData& bulkData = *mesh.getBulkData();
375
376 // loop over each entitiy extracting elements and local entity ID that
377 // are containted in specified block.
378 std::size_t entityIndex =-1;
379 std::vector<stk::mesh::Entity>::const_iterator entityItr;
380 for(entityItr=entities.begin();entityItr!=entities.end();++entityItr) {
381 stk::mesh::Entity entity = *entityItr;
382 entityIndex += 1;
383
384 const size_t num_rels = bulkData.num_elements(entity);
385 stk::mesh::Entity const* element_rels = bulkData.begin_elements(entity);
386 stk::mesh::ConnectivityOrdinal const* ordinals = bulkData.begin_element_ordinals(entity);
387 for(std::size_t e=0; e<num_rels; ++e) {
388 stk::mesh::Entity element = element_rels[e];
389 std::size_t entityId = ordinals[e];
390
391 // is this element in requested block
392 stk::mesh::Bucket const& bucket = bulkData.bucket(element);
393 bool inBlock = bucket.member(*blockPart);
394 bool onProc = bucket.member(*universalPart);
395 if(inBlock && onProc) {
396 // add element and Side ID to output vectors
397 elements.push_back(element);
398 localEntityIds.push_back(entityId);
399 } else if(!inBlock && (num_rels == 1)) {
400 // add index of side whose neighbor element in blockPart does not belong
401 // to the current processor
402 missingElementIndices.push_back(entityIndex);
403 }
404 }
405 }
406}
407
409 const std::string & blockId,
410 const std::vector<stk::mesh::Entity> & sides,
411 std::vector<std::size_t> & localSubcellDim,
412 std::vector<std::size_t> & localSubcellIds,
413 std::vector<stk::mesh::Entity> & elements)
414{
415 // This is the alogrithm, for computing the side element
416 // cascade. The requirements are that for a particular set of sides
417 // we compute all elements and subcells where they touch the side. Note
418 // that elements can be and will be repeated within this list.
419
420 std::vector<std::vector<stk::mesh::Entity> > subcells;
421 getSubcellEntities(mesh,sides,subcells);
422 subcells.push_back(sides);
423
424 // subcells now contains a unique list of faces, edges and nodes that
425 // intersect with the sides
426
427 for(std::size_t d=0;d<subcells.size();d++) {
428 std::vector<std::size_t> subcellIds;
429 std::vector<stk::mesh::Entity> subcellElements;
430
431 // find elements connected to the subcells and their local subcell information
432 getSubcellElements(mesh,blockId,subcells[d],subcellIds,subcellElements);
433
434 // sanity check
435 TEUCHOS_ASSERT(subcellIds.size()==subcellElements.size());
436
437 // concatenate with found elements
438 localSubcellDim.insert(localSubcellDim.end(),subcellElements.size(),d);
439 localSubcellIds.insert(localSubcellIds.end(),subcellIds.begin(),subcellIds.end());
440 elements.insert(elements.end(),subcellElements.begin(),subcellElements.end());
441 }
442}
443
445 const std::string & blockId,
446 const std::vector<stk::mesh::Entity> & sides,
447 std::vector<std::size_t> & localSideIds,
448 std::vector<stk::mesh::Entity> & elements)
449{
450 getSubcellElements(mesh,blockId,sides,localSideIds,elements);
451}
452
454 const std::string & blockId_a,
455 const std::string & blockId_b,
456 const std::vector<stk::mesh::Entity> & sides,
457 std::vector<std::size_t> & localSideIds_a,
458 std::vector<stk::mesh::Entity> & elements_a,
459 std::vector<std::size_t> & localSideIds_b,
460 std::vector<stk::mesh::Entity> & elements_b)
461{
462 // for verifying that an element is in specified block
463 stk::mesh::Part * blockPart_a = mesh.getElementBlockPart(blockId_a);
464 stk::mesh::Part * blockPart_b = mesh.getElementBlockPart(blockId_b);
465 stk::mesh::Part * ownedPart = mesh.getOwnedPart();
466 stk::mesh::Part * universalPart = &mesh.getMetaData()->universal_part();
467 stk::mesh::BulkData& bulkData = *mesh.getBulkData();
468
469 // loop over each entitiy extracting elements and local entity ID that
470 // are containted in specified block.
471 std::vector<stk::mesh::Entity>::const_iterator sidesItr;
472 for(sidesItr=sides.begin();sidesItr!=sides.end();++sidesItr) {
473 stk::mesh::Entity side = *sidesItr;
474
475 // these are used below the loop to insert into the appropriate vectors
476 stk::mesh::Entity element_a = stk::mesh::Entity(), element_b = stk::mesh::Entity();
477 std::size_t entityId_a=0, entityId_b=0;
478
479 const size_t num_rels = bulkData.num_elements(side);
480 stk::mesh::Entity const* element_rels = bulkData.begin_elements(side);
481 stk::mesh::ConnectivityOrdinal const* ordinals = bulkData.begin_element_ordinals(side);
482 for(std::size_t e=0; e<num_rels; ++e) {
483 stk::mesh::Entity element = element_rels[e];
484 std::size_t entityId = ordinals[e];
485
486 // is this element in requested block
487 stk::mesh::Bucket const& bucket = bulkData.bucket(element);
488 bool inBlock_a = bucket.member(*blockPart_a);
489 bool inBlock_b = bucket.member(*blockPart_b);
490 bool onProc = bucket.member(*ownedPart);
491 bool unProc = bucket.member(*universalPart);
492
493 if(inBlock_a && onProc) {
494 TEUCHOS_ASSERT(element_a==stk::mesh::Entity()); // sanity check
495 element_a = element;
496 entityId_a = entityId;
497 }
498 if(inBlock_b && unProc) {
499 TEUCHOS_ASSERT(element_b==stk::mesh::Entity()); // sanity check
500 element_b = element;
501 entityId_b = entityId;
502 }
503 }
504
505 if(element_a!=stk::mesh::Entity() && element_b!=stk::mesh::Entity()) { // add element and Side ID to output vectors
506 elements_a.push_back(element_a);
507 localSideIds_a.push_back(entityId_a);
508
509 // add element and Side ID to output vectors
510 elements_b.push_back(element_b);
511 localSideIds_b.push_back(entityId_b);
512 }
513 }
514}
515
517 const std::string & blockId,
518 const std::vector<stk::mesh::Entity> & nodes,
519 std::vector<std::size_t> & localNodeIds,
520 std::vector<stk::mesh::Entity> & elements)
521{
522 getSubcellElements(mesh,blockId,nodes,localNodeIds,elements);
523}
524
526 const std::vector<stk::mesh::Entity> & entities,
527 std::vector<std::vector<stk::mesh::Entity> > & subcells)
528{
529 // exit if there is no work to do
530 if(entities.size()==0) {
531 subcells.clear();
532 return;
533 }
534
535 stk::mesh::BulkData& bulkData = *mesh.getBulkData();
536 stk::mesh::EntityRank master_rank = bulkData.entity_rank(entities[0]);
537
538 std::vector<std::set<stk::mesh::Entity> > subcells_set(master_rank);
539
540 // loop over each entitiy extracting elements and local entity ID that
541 // are containted in specified block.
542 std::vector<stk::mesh::Entity>::const_iterator entityItr;
543 for(entityItr=entities.begin();entityItr!=entities.end();++entityItr) {
544 stk::mesh::Entity entity = *entityItr;
545
546 // sanity check, enforcing that there is only one rank
547 TEUCHOS_ASSERT(bulkData.entity_rank(entity)==master_rank);
548
549 for(int i=0; i<static_cast<int>(master_rank); i++) {
550 stk::mesh::EntityRank const to_rank = static_cast<stk::mesh::EntityRank>(i);
551 const size_t num_rels = bulkData.num_connectivity(entity, to_rank);
552 stk::mesh::Entity const* relations = bulkData.begin(entity, to_rank);
553
554 // for each relation insert the appropriate entity (into the set
555 // which gurantees uniqueness
556 for(std::size_t e=0; e<num_rels; ++e) {
557 stk::mesh::Entity subcell = relations[e];
558
559 subcells_set[i].insert(subcell);
560 }
561 }
562 }
563
564 // copy unique entities into vector
565 subcells.clear();
566 subcells.resize(subcells_set.size());
567 for(std::size_t i=0;i<subcells_set.size();i++)
568 subcells[i].assign(subcells_set[i].begin(),subcells_set[i].end());
569}
570
571}
572
573}
int baseCellDimension() const
Dimension of the base cell. NOT the dimension of the local side, even if the side() method returns tr...
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
void getElementNodes(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
std::size_t elementLocalId(stk::mesh::Entity elmt) const
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
void getMySides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
void getElementBlockNames(std::vector< std::string > &names) const
void getSidesetNames(std::vector< std::string > &name) const
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
stk::mesh::EntityRank getSideRank() const
stk::mesh::Part * getOwnedPart() const
Get a pointer to the locally owned part.
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
stk::mesh::Part * getSideset(const std::string &name) const
void getSideElements(const panzer_stk::STK_Interface &mesh, const std::string &blockId, const std::vector< stk::mesh::Entity > &sides, std::vector< std::size_t > &localSideIds, std::vector< stk::mesh::Entity > &elements)
void getUniversalSubcellElements(const panzer_stk::STK_Interface &mesh, const std::string &blockId, const std::vector< stk::mesh::Entity > &entities, std::vector< std::size_t > &localEntityIds, std::vector< stk::mesh::Entity > &elements, std::vector< std::size_t > &missingElementIndices)
void getSubcellEntities(const panzer_stk::STK_Interface &mesh, const std::vector< stk::mesh::Entity > &entities, std::vector< std::vector< stk::mesh::Entity > > &subcells)
void getNodeElements(const panzer_stk::STK_Interface &mesh, const std::string &blockId, const std::vector< stk::mesh::Entity > &nodes, std::vector< std::size_t > &localNodeIds, std::vector< stk::mesh::Entity > &elements)
void getSideElementCascade(const panzer_stk::STK_Interface &mesh, const std::string &blockId, const std::vector< stk::mesh::Entity > &sides, std::vector< std::size_t > &localSubcellDim, std::vector< std::size_t > &localSubcellIds, std::vector< stk::mesh::Entity > &elements)
void getSubcellElements(const panzer_stk::STK_Interface &mesh, const std::string &blockId, const std::vector< stk::mesh::Entity > &entities, std::vector< std::size_t > &localEntityIds, std::vector< stk::mesh::Entity > &elements)
Teuchos::RCP< std::map< unsigned, panzer::Workset > > buildBCWorksets(const panzer_stk::STK_Interface &mesh, const panzer::WorksetNeeds &needs_a, const std::string &blockid_a, const panzer::WorksetNeeds &needs_b, const std::string &blockid_b, const std::string &sideset)
Teuchos::RCP< std::vector< panzer::Workset > > buildWorksets(const panzer_stk::STK_Interface &mesh, const std::string &eBlock, const panzer::WorksetNeeds &needs)
Teuchos::RCP< std::vector< Workset > > buildWorksets(const WorksetNeeds &needs, const std::string &elementBlock, const std::vector< std::size_t > &local_cell_ids, const ArrayT &node_coordinates)
Teuchos::RCP< std::map< unsigned, Workset > > buildBCWorkset(const WorksetNeeds &needs, const std::string &elementBlock, const std::vector< std::size_t > &local_cell_ids, const std::vector< std::size_t > &local_side_ids, const ArrayT &node_coordinates, const bool populate_value_arrays=true)