12 #include "Panzer_Workset_Builder.hpp"
13 #include "Teuchos_Assert.hpp"
15 #include <stk_mesh/base/Selector.hpp>
16 #include <stk_mesh/base/GetEntities.hpp>
18 namespace panzer_stk {
22 const std::string & eBlock,
25 using namespace workset_utils;
27 std::vector<std::string> element_blocks;
29 std::vector<std::size_t> local_cell_ids;
30 Kokkos::DynRankView<double,PHX::Device> cell_node_coordinates;
32 getIdsAndNodes(mesh, eBlock, local_cell_ids, cell_node_coordinates);
44 const std::string & sideset,
45 const std::string & eBlock,
48 using namespace workset_utils;
51 std::vector<stk::mesh::Entity> sideEntities;
63 std::vector<std::string> sideSets;
67 ss << e.what() <<
"\nChoose one of:\n";
68 for(std::size_t i=0;i<sideSets.size();i++)
69 ss <<
"\"" << sideSets[i] <<
"\"\n";
75 std::vector<std::string> elementBlocks;
79 ss << e.what() <<
"\nChoose one of:\n";
80 for(std::size_t i=0;i<elementBlocks.size();i++)
81 ss <<
"\"" << elementBlocks[i] <<
"\"\n";
85 catch(std::logic_error & e) {
87 ss << e.what() <<
"\nUnrecognized logic error.\n";
92 std::vector<stk::mesh::Entity> elements;
93 std::map<std::pair<unsigned,unsigned>,std::vector<std::size_t> > local_cell_ids;
96 std::vector<std::size_t> local_side_ids;
98 sideEntities,local_side_ids,elements);
101 for(std::size_t elm=0;elm<elements.size();++elm) {
102 stk::mesh::Entity element = elements[elm];
104 local_cell_ids[std::make_pair(subcell_dim,local_side_ids[elm])].push_back(mesh.
elementLocalId(element));
108 std::vector<std::size_t> local_subcell_ids, subcell_dim;
110 sideEntities,subcell_dim,local_subcell_ids,elements);
113 for(std::size_t elm=0;elm<elements.size();++elm) {
114 stk::mesh::Entity element = elements[elm];
116 local_cell_ids[std::make_pair(subcell_dim[elm],local_subcell_ids[elm])].push_back(mesh.
elementLocalId(element));
124 if(elements.size()!=0) {
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) {
134 if(itr->second.size()==0)
137 Kokkos::DynRankView<double,PHX::Device> nodes;
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;
150 worksets->insert(worksets->end(),current->begin(),current->end());
157 return Teuchos::rcp(
new std::vector<panzer::Workset>());
163 const std::string & blockid_a,
165 const std::string & blockid_b,
166 const std::string & sideset)
168 using namespace workset_utils;
171 std::vector<stk::mesh::Entity> sideEntities;
181 stk::mesh::Part * sidePart = mesh.
getSideset(sideset);
183 "Unknown side set \"" << sideset <<
"\"");
185 stk::mesh::Selector side = *sidePart;
192 std::stringstream ss;
193 std::vector<std::string> elementBlocks;
197 ss << e.what() <<
"\nChoose one of:\n";
198 for(std::size_t i=0;i<elementBlocks.size();i++)
199 ss <<
"\"" << elementBlocks[i] <<
"\"\n";
203 catch(std::logic_error & e) {
204 std::stringstream ss;
205 ss << e.what() <<
"\nUnrecognized logic error.\n";
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;
216 local_side_ids_a,elements_a,
217 local_side_ids_b,elements_b);
220 "For a DG type boundary, the number of elements on the \"left\" and \"right\" is not the same.");
226 if(elements_a.size()==0)
227 return Teuchos::rcp(
new std::map<unsigned,panzer::Workset>);
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];
240 Kokkos::DynRankView<double,PHX::Device> node_coordinates_a, node_coordinates_b;
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);
252 const std::string & eblockID,
253 const std::string & sidesetID)
255 using namespace workset_utils;
258 std::vector<stk::mesh::Entity> sideEntities;
266 std::stringstream ss;
267 std::vector<std::string> sideSets;
271 ss << e.what() <<
"\nChoose one of:\n";
272 for(std::size_t i=0;i<sideSets.size();i++)
273 ss <<
"\"" << sideSets[i] <<
"\"\n";
278 std::stringstream ss;
279 std::vector<std::string> elementBlocks;
283 ss << e.what() <<
"\nChoose one of:\n";
284 for(std::size_t i=0;i<elementBlocks.size();i++)
285 ss <<
"\"" << elementBlocks[i] <<
"\"\n";
289 catch(std::logic_error & e) {
290 std::stringstream ss;
291 ss << e.what() <<
"\nUnrecognized logic error.\n";
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;
302 for(std::size_t elm=0;elm<elements.size();++elm) {
303 stk::mesh::Entity element = elements[elm];
312 if(elements.size()!=0) {
316 Kokkos::DynRankView<double,PHX::Device> nodes;
322 return Teuchos::null;
325 namespace workset_utils {
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)
336 stk::mesh::BulkData& bulkData = *mesh.
getBulkData();
340 std::vector<stk::mesh::Entity>::const_iterator entityItr;
341 for(entityItr=entities.begin();entityItr!=entities.end();++entityItr) {
342 stk::mesh::Entity entity = *entityItr;
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];
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) {
357 elements.push_back(element);
358 localEntityIds.push_back(entityId);
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)
373 stk::mesh::Part * universalPart = &mesh.
getMetaData()->universal_part();
374 stk::mesh::BulkData& bulkData = *mesh.
getBulkData();
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;
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];
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) {
397 elements.push_back(element);
398 localEntityIds.push_back(entityId);
399 }
else if(!inBlock && (num_rels == 1)) {
402 missingElementIndices.push_back(entityIndex);
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)
420 std::vector<std::vector<stk::mesh::Entity> > subcells;
422 subcells.push_back(sides);
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;
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());
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)
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)
466 stk::mesh::Part * universalPart = &mesh.
getMetaData()->universal_part();
467 stk::mesh::BulkData& bulkData = *mesh.
getBulkData();
471 std::vector<stk::mesh::Entity>::const_iterator sidesItr;
472 for(sidesItr=sides.begin();sidesItr!=sides.end();++sidesItr) {
473 stk::mesh::Entity side = *sidesItr;
476 stk::mesh::Entity element_a = stk::mesh::Entity(), element_b = stk::mesh::Entity();
477 std::size_t entityId_a=0, entityId_b=0;
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];
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);
493 if(inBlock_a && onProc) {
496 entityId_a = entityId;
498 if(inBlock_b && unProc) {
501 entityId_b = entityId;
505 if(element_a!=stk::mesh::Entity() && element_b!=stk::mesh::Entity()) {
506 elements_a.push_back(element_a);
507 localSideIds_a.push_back(entityId_a);
510 elements_b.push_back(element_b);
511 localSideIds_b.push_back(entityId_b);
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)
526 const std::vector<stk::mesh::Entity> & entities,
527 std::vector<std::vector<stk::mesh::Entity> > & subcells)
530 if(entities.size()==0) {
535 stk::mesh::BulkData& bulkData = *mesh.
getBulkData();
536 stk::mesh::EntityRank master_rank = bulkData.entity_rank(entities[0]);
538 std::vector<std::set<stk::mesh::Entity> > subcells_set(master_rank);
542 std::vector<stk::mesh::Entity>::const_iterator entityItr;
543 for(entityItr=entities.begin();entityItr!=entities.end();++entityItr) {
544 stk::mesh::Entity entity = *entityItr;
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);
556 for(std::size_t e=0; e<num_rels; ++e) {
557 stk::mesh::Entity subcell = relations[e];
559 subcells_set[i].insert(subcell);
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());
void getSidesetNames(std::vector< std::string > &name) const
int baseCellDimension() const
Dimension of the base cell. NOT the dimension of the local side, even if the side() method returns tr...
void getElementBlockNames(std::vector< std::string > &names) const
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
void getSubcellEntities(const panzer_stk::STK_Interface &mesh, const std::vector< stk::mesh::Entity > &entities, std::vector< std::vector< stk::mesh::Entity > > &subcells)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void getElementNodes(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
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::vector< panzer::Workset > > buildWorksets(const panzer_stk::STK_Interface &mesh, const std::string &eBlock, const panzer::WorksetNeeds &needs)
void getIdsAndNodes(const panzer_stk::STK_Interface &mesh, std::string blockId, std::vector< std::size_t > &localIds, ArrayT &nodes)
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)
stk::mesh::Part * getSideset(const std::string &name) const
void getMySides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
std::size_t elementLocalId(stk::mesh::Entity elmt) 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)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
stk::mesh::EntityRank getSideRank() const
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)
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
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)
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
#define TEUCHOS_ASSERT(assertion_test)
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 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)
stk::mesh::Part * getOwnedPart() const
Get a pointer to the locally owned part.
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
#define TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(throw_exception_test, Exception, msg)