62#include "Phalanx_KokkosDeviceTypes.hpp"
64#include "Teuchos_Assert.hpp"
65#include "Teuchos_OrdinalTraits.hpp"
67#include "Tpetra_Import.hpp"
72#include <unordered_set>
85Kokkos::DynRankView<double,PHX::Device>
86buildGhostedVertices(
const Tpetra::Import<int,panzer::GlobalOrdinal,panzer::TpetraNodeType> & importer,
87 Kokkos::DynRankView<const double,PHX::Device> owned_vertices)
92 typedef Tpetra::MultiVector<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> mvec_type;
93 typedef typename mvec_type::dual_view_type dual_view_type;
95 size_t owned_cell_cnt = importer.getSourceMap()->getLocalNumElements();
96 size_t ghstd_cell_cnt = importer.getTargetMap()->getLocalNumElements();
97 int vertices_per_cell = owned_vertices.extent(1);
98 int space_dim = owned_vertices.extent(2);
100 TEUCHOS_ASSERT(owned_vertices.extent(0)==owned_cell_cnt);
103 RCP<mvec_type> owned_vertices_mv = rcp(
new mvec_type(importer.getSourceMap(),vertices_per_cell*space_dim));
104 RCP<mvec_type> ghstd_vertices_mv = rcp(
new mvec_type(importer.getTargetMap(),vertices_per_cell*space_dim));
107 auto owned_vertices_view = owned_vertices_mv->getLocalViewDevice(Tpetra::Access::OverwriteAll);
108 Kokkos::parallel_for(owned_cell_cnt, KOKKOS_LAMBDA (
size_t i) {
110 for(
int j=0;j<vertices_per_cell;j++)
111 for(
int k=0;k<space_dim;k++,l++)
112 owned_vertices_view(i,l) = owned_vertices(i,j,k);
117 ghstd_vertices_mv->doImport(*owned_vertices_mv,importer,Tpetra::INSERT);
120 Kokkos::DynRankView<double,PHX::Device> ghstd_vertices(
"ghstd_vertices",ghstd_cell_cnt,vertices_per_cell,space_dim);
122 auto ghstd_vertices_view = ghstd_vertices_mv->getLocalViewDevice(Tpetra::Access::ReadOnly);
123 Kokkos::parallel_for(ghstd_cell_cnt, KOKKOS_LAMBDA (
size_t i) {
125 for(
int j=0;j<vertices_per_cell;j++)
126 for(
int k=0;k<space_dim;k++,l++)
127 ghstd_vertices(i,j,k) = ghstd_vertices_view(i,l);
132 return ghstd_vertices;
139 const std::string & element_block_name,
150 const shards::CellTopology & topology = *(mesh.
getCellTopology(element_block_name));
151 Teuchos::RCP<panzer::FieldPattern> cell_pattern;
152 if(topology.getDimension() == 1){
154 }
else if(topology.getDimension() == 2){
156 }
else if(topology.getDimension() == 3){
161 PANZER_FUNC_TIME_MONITOR(
"Build connectivity");
166 std::vector<panzer::LocalOrdinal> owned_block_cells;
167 auto local_cells_h = Kokkos::create_mirror_view(mesh_info.
local_cells);
168 Kokkos::deep_copy(local_cells_h, mesh_info.
local_cells);
169 for(
int parent_owned_cell=0;parent_owned_cell<num_parent_owned_cells;++parent_owned_cell){
170 const panzer::LocalOrdinal local_cell = local_cells_h(parent_owned_cell);
171 const bool is_in_block = conn.
getBlockId(local_cell) == element_block_name;
174 owned_block_cells.push_back(parent_owned_cell);
179 if ( owned_block_cells.size() == 0 )
185 PANZER_FUNC_TIME_MONITOR(
"panzer::partitioning_utilities::setupSubLocalMeshInfo");
195 const std::string & element_block_name,
196 const std::string & sideset_name,
201 using LocalOrdinal = panzer::LocalOrdinal;
202 using ParentOrdinal = int;
203 using SubcellOrdinal = short;
209 std::vector<stk::mesh::Entity> side_entities;
215 mesh.
getAllSides(sideset_name, element_block_name, side_entities);
217 }
catch(STK_Interface::SidesetException & e) {
218 std::stringstream ss;
219 std::vector<std::string> sideset_names;
223 ss << e.what() <<
"\nChoose existing sideset:\n";
224 for(
const auto & optional_sideset_name : sideset_names){
225 ss <<
"\t\"" << optional_sideset_name <<
"\"\n";
228 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true,std::logic_error,ss.str());
230 }
catch(STK_Interface::ElementBlockException & e) {
231 std::stringstream ss;
232 std::vector<std::string> element_block_names;
236 ss << e.what() <<
"\nChoose existing element block:\n";
237 for(
const auto & optional_element_block_name : element_block_names){
238 ss <<
"\t\"" << optional_element_block_name <<
"\"\n";
241 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true,std::logic_error,ss.str());
243 }
catch(std::logic_error & e) {
244 std::stringstream ss;
245 ss << e.what() <<
"\nUnrecognized logic error.\n";
247 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(
true,std::logic_error,ss.str());
252 std::map<ParentOrdinal,std::vector<SubcellOrdinal> > owned_parent_cell_to_subcell_indexes;
256 const size_t face_subcell_dimension =
static_cast<size_t>(mesh.
getCellTopology(element_block_name)->getDimension()-1);
265 std::vector<stk::mesh::Entity> elements;
266 std::vector<size_t> subcell_indexes;
267 std::vector<size_t> subcell_dimensions;
269 const size_t num_elements = subcell_dimensions.size();
272 std::unordered_map<LocalOrdinal,ParentOrdinal> local_to_parent_lookup;
273 auto local_cells_h = Kokkos::create_mirror_view(mesh_info.
local_cells);
274 Kokkos::deep_copy(local_cells_h, mesh_info.
local_cells);
275 for(ParentOrdinal parent_index=0; parent_index<static_cast<ParentOrdinal>(mesh_info.
local_cells.extent(0)); ++parent_index)
276 local_to_parent_lookup[local_cells_h(parent_index)] = parent_index;
280 for(
size_t element_index=0; element_index<num_elements; ++element_index) {
281 const size_t subcell_dimension = subcell_dimensions[element_index];
284 if(subcell_dimension == face_subcell_dimension){
285 const SubcellOrdinal subcell_index =
static_cast<SubcellOrdinal
>(subcell_indexes[element_index]);
286 const LocalOrdinal element_local_index =
static_cast<LocalOrdinal
>(mesh.
elementLocalId(elements[element_index]));
289 const auto itr = local_to_parent_lookup.find(element_local_index);
290 TEUCHOS_ASSERT(itr!= local_to_parent_lookup.end());
291 const ParentOrdinal element_parent_index = itr->second;
293 owned_parent_cell_to_subcell_indexes[element_parent_index].push_back(subcell_index);
300 const panzer::LocalOrdinal num_owned_cells = owned_parent_cell_to_subcell_indexes.size();
309 face_t(
const ParentOrdinal c0,
310 const ParentOrdinal c1,
311 const SubcellOrdinal sc0,
312 const SubcellOrdinal sc1)
319 ParentOrdinal cell_0;
320 ParentOrdinal cell_1;
321 SubcellOrdinal subcell_index_0;
322 SubcellOrdinal subcell_index_1;
327 std::unordered_set<panzer::LocalOrdinal> owned_parent_cells_set, ghstd_parent_cells_set, virtual_parent_cells_set;
328 std::vector<face_t> faces;
330 auto cell_to_faces_h = Kokkos::create_mirror_view(mesh_info.
cell_to_faces);
331 auto face_to_cells_h = Kokkos::create_mirror_view(mesh_info.
face_to_cells);
332 auto face_to_lidx_h = Kokkos::create_mirror_view(mesh_info.
face_to_lidx);
335 Kokkos::deep_copy(face_to_lidx_h, mesh_info.
face_to_lidx);
337 for(
const auto & local_cell_index_pair : owned_parent_cell_to_subcell_indexes){
339 const ParentOrdinal parent_cell = local_cell_index_pair.first;
340 const auto & subcell_indexes = local_cell_index_pair.second;
342 owned_parent_cells_set.insert(parent_cell);
344 for(
const SubcellOrdinal & subcell_index : subcell_indexes){
346 const LocalOrdinal face = cell_to_faces_h(parent_cell, subcell_index);
347 const LocalOrdinal face_other_side = (face_to_cells_h(face,0) == parent_cell) ? 1 : 0;
349 TEUCHOS_ASSERT(subcell_index == face_to_lidx_h(face, 1-face_other_side));
351 const ParentOrdinal other_side_cell = face_to_cells_h(face, face_other_side);
352 const SubcellOrdinal other_side_subcell_index = face_to_lidx_h(face, face_other_side);
354 faces.push_back(face_t(parent_cell, other_side_cell, subcell_index, other_side_subcell_index));
356 if(other_side_cell >= parent_virtual_cell_offset){
357 virtual_parent_cells_set.insert(other_side_cell);
359 ghstd_parent_cells_set.insert(other_side_cell);
365 std::vector<ParentOrdinal> all_cells;
366 all_cells.insert(all_cells.end(),owned_parent_cells_set.begin(),owned_parent_cells_set.end());
367 all_cells.insert(all_cells.end(),ghstd_parent_cells_set.begin(),ghstd_parent_cells_set.end());
368 all_cells.insert(all_cells.end(),virtual_parent_cells_set.begin(),virtual_parent_cells_set.end());
374 const LocalOrdinal num_total_cells = num_real_cells + sideset_info.
num_virtual_cells;
375 const LocalOrdinal num_vertices_per_cell = mesh_info.
cell_vertices.extent(1);
376 const LocalOrdinal num_dims = mesh_info.
cell_vertices.extent(2);
380 sideset_info.
global_cells = PHX::View<panzer::GlobalOrdinal*>(
"global_cells", num_total_cells);
381 sideset_info.
local_cells = PHX::View<LocalOrdinal*>(
"local_cells", num_total_cells);
382 sideset_info.
cell_vertices = PHX::View<double***>(
"cell_vertices", num_total_cells, num_vertices_per_cell, num_dims);
385 typename PHX::View<ParentOrdinal*>::HostMirror all_cells_h(
"all_cells_h", num_total_cells);
386 PHX::View<ParentOrdinal*> all_cells_d(
"all_cells_d", num_total_cells);
387 for(LocalOrdinal i=0; i<num_total_cells; ++i)
388 all_cells_h(i) = all_cells[i];
389 Kokkos::deep_copy(all_cells_d, all_cells_h);
390 Kokkos::parallel_for(num_total_cells, KOKKOS_LAMBDA (LocalOrdinal i) {
391 const ParentOrdinal parent_cell = all_cells_d(i);
394 for(LocalOrdinal j=0; j<num_vertices_per_cell; ++j)
395 for(LocalOrdinal k=0; k<num_dims; ++k)
402 const LocalOrdinal num_faces = faces.size();
403 const LocalOrdinal num_faces_per_cell = mesh_info.
cell_to_faces.extent(1);
405 sideset_info.
face_to_cells = PHX::View<LocalOrdinal*[2]>(
"face_to_cells", num_faces);
406 sideset_info.
face_to_lidx = PHX::View<LocalOrdinal*[2]>(
"face_to_lidx", num_faces);
407 sideset_info.
cell_to_faces = PHX::View<LocalOrdinal**>(
"cell_to_faces", num_total_cells, num_faces_per_cell);
408 auto cell_to_faces_h = Kokkos::create_mirror_view(sideset_info.
cell_to_faces);
409 auto face_to_cells_h = Kokkos::create_mirror_view(sideset_info.
face_to_cells);
410 auto face_to_lidx_h = Kokkos::create_mirror_view(sideset_info.
face_to_lidx);
413 Kokkos::deep_copy(cell_to_faces_h, -1);
416 std::unordered_map<ParentOrdinal,ParentOrdinal> parent_to_sideset_lookup;
417 for(
unsigned int i=0; i<all_cells.size(); ++i)
418 parent_to_sideset_lookup[all_cells[i]] = i;
420 for(
int face_index=0;face_index<num_faces;++face_index){
421 const face_t & face = faces[face_index];
422 const ParentOrdinal & parent_cell_0 = face.cell_0;
423 const ParentOrdinal & parent_cell_1 = face.cell_1;
426 const auto itr_0 = parent_to_sideset_lookup.find(parent_cell_0);
427 TEUCHOS_ASSERT(itr_0 != parent_to_sideset_lookup.end());
428 const ParentOrdinal sideset_cell_0 = itr_0->second;
430 const auto itr_1 = parent_to_sideset_lookup.find(parent_cell_1);
431 TEUCHOS_ASSERT(itr_1 != parent_to_sideset_lookup.end());
432 const ParentOrdinal sideset_cell_1 = itr_1->second;
437 face_to_cells_h(face_index,0) = sideset_cell_0;
438 face_to_cells_h(face_index,1) = sideset_cell_1;
440 face_to_lidx_h(face_index,0) = face.subcell_index_0;
441 face_to_lidx_h(face_index,1) = face.subcell_index_1;
443 cell_to_faces_h(sideset_cell_0,face.subcell_index_0) = face_index;
444 cell_to_faces_h(sideset_cell_1,face.subcell_index_1) = face_index;
447 Kokkos::deep_copy(sideset_info.
cell_to_faces, cell_to_faces_h);
448 Kokkos::deep_copy(sideset_info.
face_to_cells, face_to_cells_h);
449 Kokkos::deep_copy(sideset_info.
face_to_lidx, face_to_lidx_h );
455Teuchos::RCP<panzer::LocalMeshInfo>
462 typedef Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType> map_type;
463 typedef Tpetra::Import<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType> import_type;
468 auto & mesh_info = *mesh_info_rcp;
474 TEUCHOS_ASSERT(
typeid(panzer::LocalOrdinal) ==
typeid(
int));
476 Teuchos::RCP<const Teuchos::Comm<int> > comm = mesh.
getComm();
478 TEUCHOS_FUNC_TIME_MONITOR_DIFF(
"panzer_stk::generateLocalMeshInfo",GenerateLocalMeshInfo);
481 RCP<const panzer_stk::STK_Interface> mesh_rcp = Teuchos::rcpFromRef(mesh);
488 PHX::View<panzer::GlobalOrdinal*> owned_cells, ghost_cells, virtual_cells;
494 RCP<map_type> owned_cell_map = rcp(
new map_type(Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(),owned_cells,0,comm));
495 RCP<map_type> ghstd_cell_map = rcp(
new map_type(Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(),ghost_cells,0,comm));
498 RCP<import_type> cellimport_own2ghst = rcp(
new import_type(owned_cell_map,ghstd_cell_map));
504 std::vector<std::string> element_block_names;
507 const shards::CellTopology & cell_topology = *(mesh.
getCellTopology(element_block_names[0]));
509 const int space_dim = cell_topology.getDimension();
510 const int vertices_per_cell = cell_topology.getVertexCount();
511 const int faces_per_cell = cell_topology.getSubcellCount(space_dim-1);
514 Kokkos::DynRankView<double,PHX::Device> owned_vertices(
"owned_vertices",owned_cells.extent(0),vertices_per_cell,space_dim);
516 std::vector<std::size_t> localCells(owned_cells.extent(0),Teuchos::OrdinalTraits<std::size_t>::invalid());
517 for(
size_t i=0;i<localCells.size();i++)
523 Kokkos::DynRankView<double,PHX::Device> ghstd_vertices = buildGhostedVertices(*cellimport_own2ghst,owned_vertices);
528 auto owned_cells_h = Kokkos::create_mirror_view(owned_cells);
529 auto ghost_cells_h = Kokkos::create_mirror_view(ghost_cells);
530 Kokkos::deep_copy(owned_cells_h, owned_cells);
531 Kokkos::deep_copy(ghost_cells_h, ghost_cells);
532 std::unordered_map<panzer::GlobalOrdinal,int> global_to_local;
533 global_to_local[-1] = -1;
534 for(
size_t i=0;i<owned_cells.extent(0);i++)
535 global_to_local[owned_cells_h(i)] = i;
536 for(
size_t i=0;i<ghost_cells.extent(0);i++)
537 global_to_local[ghost_cells_h(i)] = i+Teuchos::as<int>(owned_cells.extent(0));
541 faceToElement->initialize(conn, comm);
542 auto elems_by_face = faceToElement->getFaceToElementsMap();
543 auto face_to_lidx = faceToElement->getFaceToCellLocalIdxMap();
555 const panzer::LocalOrdinal num_owned_cells = owned_cells.extent(0);
556 const panzer::LocalOrdinal num_ghstd_cells = ghost_cells.extent(0);
557 const panzer::LocalOrdinal num_virtual_cells = virtual_cells.extent(0);
560 const panzer::LocalOrdinal num_real_cells = num_owned_cells + num_ghstd_cells;
561 const panzer::LocalOrdinal num_total_cells = num_real_cells + num_virtual_cells;
562 const panzer::LocalOrdinal num_total_faces = elems_by_face.extent(0);
565 PHX::View<panzer::LocalOrdinal*[2]> face_to_cells(
"face_to_cells",num_total_faces);
568 PHX::View<panzer::LocalOrdinal*[2]> face_to_localidx(
"face_to_localidx",num_total_faces);
571 PHX::View<panzer::LocalOrdinal**> cell_to_face(
"cell_to_face",num_total_cells,faces_per_cell);
574 Kokkos::deep_copy(cell_to_face,-1);
578 PANZER_FUNC_TIME_MONITOR_DIFF(
"Transer faceToElement to local",TransferFaceToElementLocal);
580 int virtual_cell_index = num_real_cells;
581 auto elems_by_face_h = Kokkos::create_mirror_view(elems_by_face);
582 auto face_to_lidx_h = Kokkos::create_mirror_view(face_to_lidx);
583 auto face_to_cells_h = Kokkos::create_mirror_view(face_to_cells);
584 auto face_to_localidx_h = Kokkos::create_mirror_view(face_to_localidx);
585 auto cell_to_face_h = Kokkos::create_mirror_view(cell_to_face);
586 Kokkos::deep_copy(elems_by_face_h, elems_by_face);
587 Kokkos::deep_copy(face_to_lidx_h, face_to_lidx);
588 for(
size_t f=0;f<elems_by_face.extent(0);f++) {
590 const panzer::GlobalOrdinal global_c0 = elems_by_face_h(f,0);
591 const panzer::GlobalOrdinal global_c1 = elems_by_face_h(f,1);
594 TEUCHOS_ASSERT(global_to_local.find(global_c0)!=global_to_local.end());
595 TEUCHOS_ASSERT(global_to_local.find(global_c1)!=global_to_local.end());
597 auto c0 = global_to_local[global_c0];
598 auto lidx0 = face_to_lidx_h(f,0);
599 auto c1 = global_to_local[global_c1];
600 auto lidx1 = face_to_lidx_h(f,1);
607 c0 = virtual_cell_index++;
614 cell_to_face_h(c0,lidx0) = f;
620 c1 = virtual_cell_index++;
627 cell_to_face_h(c1,lidx1) = f;
631 face_to_cells_h(f,0) = c0;
632 face_to_localidx_h(f,0) = lidx0;
633 face_to_cells_h(f,1) = c1;
634 face_to_localidx_h(f,1) = lidx1;
636 face_to_cells_h(f,0) = c1;
637 face_to_localidx_h(f,0) = lidx1;
638 face_to_cells_h(f,1) = c0;
639 face_to_localidx_h(f,1) = lidx0;
643 TEUCHOS_ASSERT(c0<num_real_cells or c1<num_real_cells)
646 Kokkos::deep_copy(face_to_cells, face_to_cells_h);
647 Kokkos::deep_copy(face_to_localidx, face_to_localidx_h);
648 Kokkos::deep_copy(cell_to_face, cell_to_face_h);
655 PANZER_FUNC_TIME_MONITOR_DIFF(
"Assign Indices",AssignIndices);
667 mesh_info.
global_cells = PHX::View<panzer::GlobalOrdinal*>(
"global_cell_indices",num_total_cells);
668 mesh_info.
local_cells = PHX::View<panzer::LocalOrdinal*>(
"local_cell_indices",num_total_cells);
670 Kokkos::parallel_for(num_owned_cells,KOKKOS_LAMBDA (
int i) {
675 Kokkos::parallel_for(num_ghstd_cells,KOKKOS_LAMBDA (
int i) {
676 mesh_info.
global_cells(i+num_owned_cells) = ghost_cells(i);
677 mesh_info.
local_cells(i+num_owned_cells) = i+num_owned_cells;
680 Kokkos::parallel_for(num_virtual_cells,KOKKOS_LAMBDA (
int i) {
681 mesh_info.
global_cells(i+num_real_cells) = virtual_cells(i);
682 mesh_info.
local_cells(i+num_real_cells) = i+num_real_cells;
685 mesh_info.
cell_vertices = PHX::View<double***>(
"cell_vertices",num_total_cells,vertices_per_cell,space_dim);
690 Kokkos::parallel_for(num_owned_cells,KOKKOS_LAMBDA (
int i) {
691 for(
int j=0;j<vertices_per_cell;++j){
692 for(
int k=0;k<space_dim;++k){
698 Kokkos::parallel_for(num_ghstd_cells,KOKKOS_LAMBDA (
int i) {
699 for(
int j=0;j<vertices_per_cell;++j){
700 for(
int k=0;k<space_dim;++k){
701 mesh_info.
cell_vertices(i+num_owned_cells,j,k) = ghstd_vertices(i,j,k);
709 PANZER_FUNC_TIME_MONITOR_DIFF(
"Assign geometry traits",AssignGeometryTraits);
710 Kokkos::parallel_for(num_virtual_cells,KOKKOS_LAMBDA (
int i) {
711 const panzer::LocalOrdinal virtual_cell = i+num_real_cells;
712 for(
int local_face=0; local_face<faces_per_cell; ++local_face){
713 const panzer::LocalOrdinal face = cell_to_face(virtual_cell, local_face);
715 const panzer::LocalOrdinal other_side = (face_to_cells(face, 0) == virtual_cell) ? 1 : 0;
716 const panzer::LocalOrdinal real_cell = face_to_cells(face,other_side);
717 for(
int j=0;j<vertices_per_cell;++j){
718 for(
int k=0;k<space_dim;++k){
731 std::vector<std::string> sideset_names;
734 for(
const std::string & element_block_name : element_block_names){
735 PANZER_FUNC_TIME_MONITOR_DIFF(
"Set up setupLocalMeshBlockInfo",SetupLocalMeshBlockInfo);
737 setupLocalMeshBlockInfo(mesh, conn, mesh_info, element_block_name, block_info);
743 for(
const std::string & sideset_name : sideset_names){
744 PANZER_FUNC_TIME_MONITOR_DIFF(
"Setup LocalMeshSidesetInfo",SetupLocalMeshSidesetInfo);
746 setupLocalMeshSidesetInfo(mesh, conn, mesh_info, element_block_name, sideset_name, sideset_info);
754 return mesh_info_rcp;
Pure virtual base class for supplying mesh connectivity information to the DOF Manager.
virtual std::string getBlockId(LocalOrdinal localElmtId) const =0
virtual void buildConnectivity(const FieldPattern &fp)=0
bool isInitialized() const
Has initialize been called on this mesh object?
void getElementVerticesNoResize(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
std::size_t elementLocalId(stk::mesh::Entity elmt) const
void getAllSides(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< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
void setupSubLocalMeshInfo(const panzer::LocalMeshInfoBase &parent_info, const std::vector< panzer::LocalOrdinal > &owned_parent_cells, panzer::LocalMeshInfoBase &sub_info)
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)
Teuchos::RCP< panzer::LocalMeshInfo > generateLocalMeshInfo(const panzer_stk::STK_Interface &mesh)
Create a structure containing information about the local portion of a given element block.
void fillLocalCellIDs(const Teuchos::RCP< const Teuchos::Comm< int > > &comm, panzer::ConnManager &conn, PHX::View< panzer::GlobalOrdinal * > &owned_cells, PHX::View< panzer::GlobalOrdinal * > &ghost_cells, PHX::View< panzer::GlobalOrdinal * > &virtual_cells)
Get the owned, ghost and virtual global cell ids for this process.
std::string element_block_name
Teuchos::RCP< const shards::CellTopology > cell_topology
panzer::LocalOrdinal num_owned_cells
PHX::View< panzer::LocalOrdinal *[2]> face_to_lidx
PHX::View< panzer::GlobalOrdinal * > global_cells
panzer::LocalOrdinal num_ghstd_cells
PHX::View< double *** > cell_vertices
panzer::LocalOrdinal num_virtual_cells
PHX::View< panzer::LocalOrdinal *[2]> face_to_cells
PHX::View< panzer::LocalOrdinal ** > cell_to_faces
PHX::View< panzer::LocalOrdinal * > local_cells
std::map< std::string, std::map< std::string, LocalMeshSidesetInfo > > sidesets
std::map< std::string, LocalMeshBlockInfo > element_blocks
Teuchos::RCP< const shards::CellTopology > cell_topology
std::string element_block_name