Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_LocalMeshUtilities.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 "Panzer_NodeType.hpp"
16
17#include "Panzer_HashUtils.hpp"
21
27
29
30#include "Phalanx_KokkosDeviceTypes.hpp"
31
32#include "Teuchos_Assert.hpp"
33#include "Teuchos_OrdinalTraits.hpp"
34
35#include "Tpetra_Import.hpp"
36
37#include <string>
38#include <map>
39#include <vector>
40#include <unordered_set>
41
42namespace panzer_stk
43{
44
45// No external access
46namespace
47{
48
52Kokkos::DynRankView<double,PHX::Device>
53buildGhostedNodes(const Tpetra::Import<int,panzer::GlobalOrdinal,panzer::TpetraNodeType> & importer,
54 Kokkos::DynRankView<const double,PHX::Device> owned_nodes)
55{
56 using Teuchos::RCP;
57 using Teuchos::rcp;
58
59 typedef Tpetra::MultiVector<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> mvec_type;
60
61 size_t owned_cell_cnt = importer.getSourceMap()->getLocalNumElements();
62 size_t ghstd_cell_cnt = importer.getTargetMap()->getLocalNumElements();
63 int nodes_per_cell = owned_nodes.extent(1);
64 int space_dim = owned_nodes.extent(2);
65
66 TEUCHOS_ASSERT(owned_nodes.extent(0)==owned_cell_cnt);
67
68 // build node multivector
69 RCP<mvec_type> owned_nodes_mv = rcp(new mvec_type(importer.getSourceMap(),nodes_per_cell*space_dim));
70 RCP<mvec_type> ghstd_nodes_mv = rcp(new mvec_type(importer.getTargetMap(),nodes_per_cell*space_dim));
71
72 {
73 auto owned_nodes_view = owned_nodes_mv->getLocalViewDevice(Tpetra::Access::OverwriteAll);
74 Kokkos::parallel_for(owned_cell_cnt, KOKKOS_LAMBDA (size_t i) {
75 int l = 0;
76 for(int j=0;j<nodes_per_cell;j++)
77 for(int k=0;k<space_dim;k++,l++)
78 owned_nodes_view(i,l) = owned_nodes(i,j,k);
79 });
80 }
81
82 // communicate ghstd nodes
83 ghstd_nodes_mv->doImport(*owned_nodes_mv,importer,Tpetra::INSERT);
84
85 // copy multivector into ghstd node structure
86 Kokkos::DynRankView<double,PHX::Device> ghstd_nodes("ghstd_nodes",ghstd_cell_cnt,nodes_per_cell,space_dim);
87 {
88 auto ghstd_nodes_view = ghstd_nodes_mv->getLocalViewDevice(Tpetra::Access::ReadOnly);
89 Kokkos::parallel_for(ghstd_cell_cnt, KOKKOS_LAMBDA (size_t i) {
90 int l = 0;
91 for(int j=0;j<nodes_per_cell;j++)
92 for(int k=0;k<space_dim;k++,l++)
93 ghstd_nodes(i,j,k) = ghstd_nodes_view(i,l);
94 } );
95 Kokkos::fence();
96 }
97
98 return ghstd_nodes;
99} // end build ghstd nodes
100void
101setupLocalMeshBlockInfo(const panzer_stk::STK_Interface & mesh,
102 panzer::ConnManager & conn,
103 const panzer::LocalMeshInfo & mesh_info,
104 const std::string & element_block_name,
105 panzer::LocalMeshBlockInfo & block_info)
106{
107
108 // This function identifies all cells in mesh_info that belong to element_block_name
109 // and creates a block_info from it.
110
111 const int num_parent_owned_cells = mesh_info.num_owned_cells;
112
113 // Make sure connectivity is setup for interfaces between cells
114 {
115 const shards::CellTopology & topology = *(mesh.getCellTopology(element_block_name));
116 Teuchos::RCP<panzer::FieldPattern> cell_pattern;
117 if(topology.getDimension() == 1){
118 cell_pattern = Teuchos::rcp(new panzer::EdgeFieldPattern(topology));
119 } else if(topology.getDimension() == 2){
120 cell_pattern = Teuchos::rcp(new panzer::FaceFieldPattern(topology));
121 } else if(topology.getDimension() == 3){
122 cell_pattern = Teuchos::rcp(new panzer::ElemFieldPattern(topology));
123 }
124
125 {
126 PANZER_FUNC_TIME_MONITOR("Build connectivity");
127 conn.buildConnectivity(*cell_pattern);
128 }
129 }
130
131 std::vector<panzer::LocalOrdinal> owned_block_cells;
132 auto local_cells_h = Kokkos::create_mirror_view(mesh_info.local_cells);
133 Kokkos::deep_copy(local_cells_h, mesh_info.local_cells);
134 for(int parent_owned_cell=0;parent_owned_cell<num_parent_owned_cells;++parent_owned_cell){
135 const panzer::LocalOrdinal local_cell = local_cells_h(parent_owned_cell);
136 const bool is_in_block = conn.getBlockId(local_cell) == element_block_name;
137
138 if(is_in_block){
139 owned_block_cells.push_back(parent_owned_cell);
140 }
141
142 }
143
144 if ( owned_block_cells.size() == 0 )
145 return;
146 block_info.num_owned_cells = owned_block_cells.size();
147 block_info.element_block_name = element_block_name;
148 block_info.cell_topology = mesh.getCellTopology(element_block_name);
149 {
150 PANZER_FUNC_TIME_MONITOR("panzer::partitioning_utilities::setupSubLocalMeshInfo");
151 panzer::partitioning_utilities::setupSubLocalMeshInfo(mesh_info, owned_block_cells, block_info);
152 }
153}
154
155
156void
157setupLocalMeshSidesetInfo(const panzer_stk::STK_Interface & mesh,
158 panzer::ConnManager& /* conn */,
159 const panzer::LocalMeshInfo & mesh_info,
160 const std::string & element_block_name,
161 const std::string & sideset_name,
162 panzer::LocalMeshSidesetInfo & sideset_info)
163{
164
165 // We use these typedefs to make the algorithm slightly more clear
166 using LocalOrdinal = panzer::LocalOrdinal;
167 using ParentOrdinal = int;
168 using SubcellOrdinal = short;
169
170 // This function identifies all cells in mesh_info that belong to element_block_name
171 // and creates a block_info from it.
172
173 // This is a list of all entities that make up the 'side'
174 std::vector<stk::mesh::Entity> side_entities;
175
176 // Grab the side entities associated with this sideset on the mesh
177 // Note: Throws exception if element block or sideset doesn't exist
178 try{
179
180 mesh.getAllSides(sideset_name, element_block_name, side_entities);
181
182 } catch(STK_Interface::SidesetException & e) {
183 std::stringstream ss;
184 std::vector<std::string> sideset_names;
185 mesh.getSidesetNames(sideset_names);
186
187 // build an error message
188 ss << e.what() << "\nChoose existing sideset:\n";
189 for(const auto & optional_sideset_name : sideset_names){
190 ss << "\t\"" << optional_sideset_name << "\"\n";
191 }
192
193 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
194
195 } catch(STK_Interface::ElementBlockException & e) {
196 std::stringstream ss;
197 std::vector<std::string> element_block_names;
198 mesh.getElementBlockNames(element_block_names);
199
200 // build an error message
201 ss << e.what() << "\nChoose existing element block:\n";
202 for(const auto & optional_element_block_name : element_block_names){
203 ss << "\t\"" << optional_element_block_name << "\"\n";
204 }
205
206 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
207
208 } catch(std::logic_error & e) {
209 std::stringstream ss;
210 ss << e.what() << "\nUnrecognized logic error.\n";
211
212 TEUCHOS_TEST_FOR_EXCEPTION_PURE_MSG(true,std::logic_error,ss.str());
213
214 }
215
216 // We now have a list of sideset entities, lets unwrap them and create the sideset_info!
217 std::map<ParentOrdinal,std::vector<SubcellOrdinal> > owned_parent_cell_to_subcell_indexes;
218 {
219
220 // This is the subcell dimension we are trying to line up on the sideset
221 const size_t face_subcell_dimension = static_cast<size_t>(mesh.getCellTopology(element_block_name)->getDimension()-1);
222
223 // List of local subcell indexes indexed by element:
224 // For example: a Tet (element) would have
225 // - 4 triangular faces (subcell_index 0-3, subcell_dimension=2)
226 // - 6 edges (subcell_index 0-5, subcell_dimension=1)
227 // - 4 nodes (subcell_index 0-3, subcell_dimension=0)
228 // Another example: a Line (element) would have
229 // - 2 nodes (subcell_index 0-1, subcell_dimension=0)
230 // The nodes coincide with the element vertices for these first order examples
231 std::vector<stk::mesh::Entity> elements;
232 std::vector<size_t> subcell_indexes;
233 std::vector<size_t> subcell_dimensions;
234 panzer_stk::workset_utils::getSideElementCascade(mesh, element_block_name, side_entities, subcell_dimensions, subcell_indexes, elements);
235 const size_t num_elements = subcell_dimensions.size();
236
237 // We need a fast lookup for maping local indexes to parent indexes
238 std::unordered_map<LocalOrdinal,ParentOrdinal> local_to_parent_lookup;
239 auto local_cells_h = Kokkos::create_mirror_view(mesh_info.local_cells);
240 Kokkos::deep_copy(local_cells_h, mesh_info.local_cells);
241 for(ParentOrdinal parent_index=0; parent_index<static_cast<ParentOrdinal>(mesh_info.local_cells.extent(0)); ++parent_index)
242 local_to_parent_lookup[local_cells_h(parent_index)] = parent_index;
243
244 // Add the subcell indexes for each element on the sideset
245 // TODO: There is a lookup call here to map from local indexes to parent indexes which slows things down. Maybe there is a way around that
246 for(size_t element_index=0; element_index<num_elements; ++element_index) {
247 const size_t subcell_dimension = subcell_dimensions[element_index];
248
249 // Add subcell to map
250 if(subcell_dimension == face_subcell_dimension){
251 const SubcellOrdinal subcell_index = static_cast<SubcellOrdinal>(subcell_indexes[element_index]);
252 const LocalOrdinal element_local_index = static_cast<LocalOrdinal>(mesh.elementLocalId(elements[element_index]));
253
254 // Look up the parent cell index using the local cell index
255 const auto itr = local_to_parent_lookup.find(element_local_index);
256 TEUCHOS_ASSERT(itr!= local_to_parent_lookup.end());
257 const ParentOrdinal element_parent_index = itr->second;
258
259 owned_parent_cell_to_subcell_indexes[element_parent_index].push_back(subcell_index);
260 }
261 }
262 }
263
264 // We now know the mapping of parent cell indexes to subcell indexes touching the sideset
265
266 const panzer::LocalOrdinal num_owned_cells = owned_parent_cell_to_subcell_indexes.size();
267
268 sideset_info.element_block_name = element_block_name;
269 sideset_info.sideset_name = sideset_name;
270 sideset_info.cell_topology = mesh.getCellTopology(element_block_name);
271
272 sideset_info.num_owned_cells = num_owned_cells;
273
274 struct face_t{
275 face_t(const ParentOrdinal c0,
276 const ParentOrdinal c1,
277 const SubcellOrdinal sc0,
278 const SubcellOrdinal sc1)
279 {
280 cell_0=c0;
281 cell_1=c1;
282 subcell_index_0=sc0;
283 subcell_index_1=sc1;
284 }
285 ParentOrdinal cell_0;
286 ParentOrdinal cell_1;
287 SubcellOrdinal subcell_index_0;
288 SubcellOrdinal subcell_index_1;
289 };
290
291
292 // Figure out how many cells on the other side of the sideset are ghost or virtual
293 std::unordered_set<panzer::LocalOrdinal> owned_parent_cells_set, ghstd_parent_cells_set, virtual_parent_cells_set;
294 std::vector<face_t> faces;
295 {
296 auto cell_to_faces_h = Kokkos::create_mirror_view(mesh_info.cell_to_faces);
297 auto face_to_cells_h = Kokkos::create_mirror_view(mesh_info.face_to_cells);
298 auto face_to_lidx_h = Kokkos::create_mirror_view(mesh_info.face_to_lidx);
299 Kokkos::deep_copy(cell_to_faces_h, mesh_info.cell_to_faces);
300 Kokkos::deep_copy(face_to_cells_h, mesh_info.face_to_cells);
301 Kokkos::deep_copy(face_to_lidx_h, mesh_info.face_to_lidx);
302 panzer::LocalOrdinal parent_virtual_cell_offset = mesh_info.num_owned_cells + mesh_info.num_ghstd_cells;
303 for(const auto & local_cell_index_pair : owned_parent_cell_to_subcell_indexes){
304
305 const ParentOrdinal parent_cell = local_cell_index_pair.first;
306 const auto & subcell_indexes = local_cell_index_pair.second;
307
308 owned_parent_cells_set.insert(parent_cell);
309
310 for(const SubcellOrdinal & subcell_index : subcell_indexes){
311
312 const LocalOrdinal face = cell_to_faces_h(parent_cell, subcell_index);
313 const LocalOrdinal face_other_side = (face_to_cells_h(face,0) == parent_cell) ? 1 : 0;
314
315 TEUCHOS_ASSERT(subcell_index == face_to_lidx_h(face, 1-face_other_side));
316
317 const ParentOrdinal other_side_cell = face_to_cells_h(face, face_other_side);
318 const SubcellOrdinal other_side_subcell_index = face_to_lidx_h(face, face_other_side);
319
320 faces.push_back(face_t(parent_cell, other_side_cell, subcell_index, other_side_subcell_index));
321
322 if(other_side_cell >= parent_virtual_cell_offset){
323 virtual_parent_cells_set.insert(other_side_cell);
324 } else {
325 ghstd_parent_cells_set.insert(other_side_cell);
326 }
327 }
328 }
329 }
330
331 std::vector<ParentOrdinal> all_cells;
332 all_cells.insert(all_cells.end(),owned_parent_cells_set.begin(),owned_parent_cells_set.end());
333 all_cells.insert(all_cells.end(),ghstd_parent_cells_set.begin(),ghstd_parent_cells_set.end());
334 all_cells.insert(all_cells.end(),virtual_parent_cells_set.begin(),virtual_parent_cells_set.end());
335
336 sideset_info.num_ghstd_cells = ghstd_parent_cells_set.size();
337 sideset_info.num_virtual_cells = virtual_parent_cells_set.size();
338
339 const LocalOrdinal num_real_cells = sideset_info.num_owned_cells + sideset_info.num_ghstd_cells;
340 const LocalOrdinal num_total_cells = num_real_cells + sideset_info.num_virtual_cells;
341 const LocalOrdinal num_nodes_per_cell = mesh_info.cell_nodes.extent(1);
342 const LocalOrdinal num_dims = mesh_info.cell_nodes.extent(2);
343
344 // Copy local indexes, global indexes, and cell nodes to sideset info
345 {
346 sideset_info.global_cells = PHX::View<panzer::GlobalOrdinal*>("global_cells", num_total_cells);
347 sideset_info.local_cells = PHX::View<LocalOrdinal*>("local_cells", num_total_cells);
348 sideset_info.cell_nodes = PHX::View<double***>("cell_nodes", num_total_cells, num_nodes_per_cell, num_dims);
349 Kokkos::deep_copy(sideset_info.cell_nodes,0.);
350
351 typename PHX::View<ParentOrdinal*>::host_mirror_type all_cells_h("all_cells_h", num_total_cells);
352 PHX::View<ParentOrdinal*> all_cells_d("all_cells_d", num_total_cells);
353 for(LocalOrdinal i=0; i<num_total_cells; ++i)
354 all_cells_h(i) = all_cells[i];
355 Kokkos::deep_copy(all_cells_d, all_cells_h);
356 Kokkos::parallel_for(num_total_cells, KOKKOS_LAMBDA (LocalOrdinal i) {
357 const ParentOrdinal parent_cell = all_cells_d(i);
358 sideset_info.local_cells(i) = mesh_info.local_cells(parent_cell);
359 sideset_info.global_cells(i) = mesh_info.global_cells(parent_cell);
360 for(LocalOrdinal j=0; j<num_nodes_per_cell; ++j)
361 for(LocalOrdinal k=0; k<num_dims; ++k)
362 sideset_info.cell_nodes(i,j,k) = mesh_info.cell_nodes(parent_cell,j,k);
363 });
364 }
365
366 // Now we have to set the connectivity for the faces.
367
368 const LocalOrdinal num_faces = faces.size();
369 const LocalOrdinal num_faces_per_cell = mesh_info.cell_to_faces.extent(1);
370
371 sideset_info.face_to_cells = PHX::View<LocalOrdinal*[2]>("face_to_cells", num_faces);
372 sideset_info.face_to_lidx = PHX::View<LocalOrdinal*[2]>("face_to_lidx", num_faces);
373 sideset_info.cell_to_faces = PHX::View<LocalOrdinal**>("cell_to_faces", num_total_cells, num_faces_per_cell);
374 auto cell_to_faces_h = Kokkos::create_mirror_view(sideset_info.cell_to_faces);
375 auto face_to_cells_h = Kokkos::create_mirror_view(sideset_info.face_to_cells);
376 auto face_to_lidx_h = Kokkos::create_mirror_view(sideset_info.face_to_lidx);
377
378 // Default the system with invalid cell index - this will be most of the entries
379 Kokkos::deep_copy(cell_to_faces_h, -1);
380
381 // Lookup for sideset cell index given parent cell index
382 std::unordered_map<ParentOrdinal,ParentOrdinal> parent_to_sideset_lookup;
383 for(unsigned int i=0; i<all_cells.size(); ++i)
384 parent_to_sideset_lookup[all_cells[i]] = i;
385
386 for(int face_index=0;face_index<num_faces;++face_index){
387 const face_t & face = faces[face_index];
388 const ParentOrdinal & parent_cell_0 = face.cell_0;
389 const ParentOrdinal & parent_cell_1 = face.cell_1;
390
391 // Convert the parent cell index into a sideset cell index
392 const auto itr_0 = parent_to_sideset_lookup.find(parent_cell_0);
393 TEUCHOS_ASSERT(itr_0 != parent_to_sideset_lookup.end());
394 const ParentOrdinal sideset_cell_0 = itr_0->second;
395
396 const auto itr_1 = parent_to_sideset_lookup.find(parent_cell_1);
397 TEUCHOS_ASSERT(itr_1 != parent_to_sideset_lookup.end());
398 const ParentOrdinal sideset_cell_1 = itr_1->second;
399
400// const ParentOrdinal sideset_cell_0 = std::distance(all_cells.begin(), std::find(all_cells.begin(), all_cells.end(), parent_cell_0));
401// const ParentOrdinal sideset_cell_1 = std::distance(all_cells.begin(), std::find(all_cells.begin()+num_owned_cells, all_cells.end(), parent_cell_1));
402
403 face_to_cells_h(face_index,0) = sideset_cell_0;
404 face_to_cells_h(face_index,1) = sideset_cell_1;
405
406 face_to_lidx_h(face_index,0) = face.subcell_index_0;
407 face_to_lidx_h(face_index,1) = face.subcell_index_1;
408
409 cell_to_faces_h(sideset_cell_0,face.subcell_index_0) = face_index;
410 cell_to_faces_h(sideset_cell_1,face.subcell_index_1) = face_index;
411
412 }
413 Kokkos::deep_copy(sideset_info.cell_to_faces, cell_to_faces_h);
414 Kokkos::deep_copy(sideset_info.face_to_cells, face_to_cells_h);
415 Kokkos::deep_copy(sideset_info.face_to_lidx, face_to_lidx_h );
416
417}
418
419} // namespace
420
421Teuchos::RCP<panzer::LocalMeshInfo>
423{
424 TEUCHOS_FUNC_TIME_MONITOR_DIFF("panzer_stk::generateLocalMeshInfo",GenerateLocalMeshInfo);
425
426 using Teuchos::RCP;
427 using Teuchos::rcp;
428
429 //typedef Tpetra::CrsMatrix<int,panzer::LocalOrdinal,panzer::GlobalOrdinal> crs_type;
430 typedef Tpetra::Map<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType> map_type;
431 typedef Tpetra::Import<panzer::LocalOrdinal,panzer::GlobalOrdinal,panzer::TpetraNodeType> import_type;
432 //typedef Tpetra::MultiVector<double,panzer::LocalOrdinal,panzer::GlobalOrdinal> mvec_type;
433 //typedef Tpetra::MultiVector<panzer::GlobalOrdinal,panzer::LocalOrdinal,panzer::GlobalOrdinal> ordmvec_type;
434
435 auto mesh_info_rcp = Teuchos::rcp(new panzer::LocalMeshInfo);
436 auto & mesh_info = *mesh_info_rcp;
437
438 // Make sure the STK interface is valid
439 TEUCHOS_ASSERT(mesh.isInitialized());
440
441 // This is required by some of the STK stuff
442 TEUCHOS_ASSERT(typeid(panzer::LocalOrdinal) == typeid(int));
443
444 Teuchos::RCP<const Teuchos::Comm<int> > comm = mesh.getComm();
445
446 // This horrible line of code is required since the connection manager only takes rcps of a mesh
447 RCP<const panzer_stk::STK_Interface> mesh_rcp = Teuchos::rcpFromRef(mesh);
448 // We're allowed to do this since the connection manager only exists in this scope... even though it is also an RCP...
449
450 // extract topology handle
451 RCP<panzer::ConnManager> conn_rcp = rcp(new panzer_stk::STKConnManager(mesh_rcp));
452 panzer::ConnManager & conn = *conn_rcp;
453
454 PHX::View<panzer::GlobalOrdinal*> owned_cells, ghost_cells, virtual_cells;
455 panzer::fillLocalCellIDs(comm, conn, owned_cells, ghost_cells, virtual_cells);
456
457 // build cell maps
459
460 RCP<map_type> owned_cell_map = rcp(new map_type(Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(),owned_cells,0,comm));
461 RCP<map_type> ghstd_cell_map = rcp(new map_type(Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid(),ghost_cells,0,comm));
462
463 // build importer: cell importer, owned to ghstd
464 RCP<import_type> cellimport_own2ghst = rcp(new import_type(owned_cell_map,ghstd_cell_map));
465
466 // read all the nodes associated with these elements, get ghstd contributions
468
469 // TODO: This all needs to be rewritten for when element blocks have different cell topologies
470 std::vector<std::string> element_block_names;
471 mesh.getElementBlockNames(element_block_names);
472
473 const shards::CellTopology & cell_topology = *(mesh.getCellTopology(element_block_names[0]));
474
475 const int space_dim = cell_topology.getDimension();
476 const int nodes_per_cell = cell_topology.getNodeCount();
477 const int faces_per_cell = cell_topology.getSubcellCount(space_dim-1);
478
479 Kokkos::DynRankView<double,PHX::Device> owned_nodes("owned_nodes",owned_cells.extent(0),nodes_per_cell,space_dim);
480 {
481 std::vector<std::size_t> localCells(owned_cells.extent(0),Teuchos::OrdinalTraits<std::size_t>::invalid());
482 for(size_t i=0;i<localCells.size();i++)
483 localCells[i] = i;
484 mesh.getElementNodesNoResize(localCells,owned_nodes);
485 }
486
487 // this builds a ghstd node array
488 Kokkos::DynRankView<double,PHX::Device> ghstd_nodes = buildGhostedNodes(*cellimport_own2ghst,owned_nodes);
489
490 // build edge to cell neighbor mapping
492
493 auto owned_cells_h = Kokkos::create_mirror_view(owned_cells);
494 auto ghost_cells_h = Kokkos::create_mirror_view(ghost_cells);
495 Kokkos::deep_copy(owned_cells_h, owned_cells);
496 Kokkos::deep_copy(ghost_cells_h, ghost_cells);
497 std::unordered_map<panzer::GlobalOrdinal,int> global_to_local;
498 global_to_local[-1] = -1; // this is the "no neighbor" flag
499 for(size_t i=0;i<owned_cells.extent(0);i++)
500 global_to_local[owned_cells_h(i)] = i;
501 for(size_t i=0;i<ghost_cells.extent(0);i++)
502 global_to_local[ghost_cells_h(i)] = i+Teuchos::as<int>(owned_cells.extent(0));
503
504 // this class comes from Mini-PIC and Matt B
505 RCP<panzer::FaceToElement<panzer::LocalOrdinal,panzer::GlobalOrdinal> > faceToElement = rcp(new panzer::FaceToElement<panzer::LocalOrdinal,panzer::GlobalOrdinal>());
506 faceToElement->initialize(conn, comm);
507 auto elems_by_face = faceToElement->getFaceToElementsMap();
508 auto face_to_lidx = faceToElement->getFaceToCellLocalIdxMap();
509
510 // We also need to consider faces that connect to cells that do not exist, but are needed for boundary conditions
511 // We dub them virtual cell since there should be no geometry associated with them, or topology really
512 // They exist only for datastorage so that they are consistent with 'real' cells from an algorithm perspective
513
514 // Each virtual face (face linked to a '-1' cell) requires a virtual cell (i.e. turn the '-1' into a virtual cell)
515 // Virtual cells are those that do not exist but are connected to an owned cell
516 // Note - in the future, ghosted cells will also need to connect to virtual cells at boundary conditions, but for the moment we will ignore this.
517
518 // Iterate over all faces and identify the faces connected to a potential virtual cell
519
520 const panzer::LocalOrdinal num_owned_cells = owned_cells.extent(0);
521 const panzer::LocalOrdinal num_ghstd_cells = ghost_cells.extent(0);
522 const panzer::LocalOrdinal num_virtual_cells = virtual_cells.extent(0);
523
524 // total cells and faces include owned, ghosted, and virtual
525 const panzer::LocalOrdinal num_real_cells = num_owned_cells + num_ghstd_cells;
526 const panzer::LocalOrdinal num_total_cells = num_real_cells + num_virtual_cells;
527 const panzer::LocalOrdinal num_total_faces = elems_by_face.extent(0);
528
529 // Lookup cells connected to a face
530 PHX::View<panzer::LocalOrdinal*[2]> face_to_cells("face_to_cells",num_total_faces);
531
532 // Lookup local face indexes given cell and left/right state (0/1)
533 PHX::View<panzer::LocalOrdinal*[2]> face_to_localidx("face_to_localidx",num_total_faces);
534
535 // Lookup face index given a cell and local face index
536 PHX::View<panzer::LocalOrdinal**> cell_to_face("cell_to_face",num_total_cells,faces_per_cell);
537
538 // Transfer information from 'faceToElement' datasets to local arrays
539 {
540 PANZER_FUNC_TIME_MONITOR_DIFF("Transer faceToElement to local",TransferFaceToElementLocal);
541
542 int virtual_cell_index = num_real_cells;
543 auto elems_by_face_h = Kokkos::create_mirror_view(elems_by_face);
544 auto face_to_lidx_h = Kokkos::create_mirror_view(face_to_lidx);
545 auto face_to_cells_h = Kokkos::create_mirror_view(face_to_cells);
546 auto face_to_localidx_h = Kokkos::create_mirror_view(face_to_localidx);
547 auto cell_to_face_h = Kokkos::create_mirror_view(cell_to_face);
548 Kokkos::deep_copy(elems_by_face_h, elems_by_face);
549 Kokkos::deep_copy(face_to_lidx_h, face_to_lidx);
550 // initialize with negative one cells that are not associated with a face
551 Kokkos::deep_copy(cell_to_face_h, -1);
552 for(size_t f=0;f<elems_by_face.extent(0);f++) {
553
554 const panzer::GlobalOrdinal global_c0 = elems_by_face_h(f,0);
555 const panzer::GlobalOrdinal global_c1 = elems_by_face_h(f,1);
556
557 // make sure that no bonus cells get in here
558 TEUCHOS_ASSERT(global_to_local.find(global_c0)!=global_to_local.end());
559 TEUCHOS_ASSERT(global_to_local.find(global_c1)!=global_to_local.end());
560
561 auto c0 = global_to_local[global_c0];
562 auto lidx0 = face_to_lidx_h(f,0);
563 auto c1 = global_to_local[global_c1];
564 auto lidx1 = face_to_lidx_h(f,1);
565
566 // Test for virtual cells
567
568 // Left cell
569 if(c0 < 0){
570 // Virtual cell - create it!
571 c0 = virtual_cell_index++;
572
573 // We need the subcell_index to line up between real and virtual cell
574 // This way the face has the same geometry... though the face normal
575 // will point in the wrong direction
576 lidx0 = lidx1;
577 }
578 cell_to_face_h(c0,lidx0) = f;
579
580
581 // Right cell
582 if(c1 < 0){
583 // Virtual cell - create it!
584 c1 = virtual_cell_index++;
585
586 // We need the subcell_index to line up between real and virtual cell
587 // This way the face has the same geometry... though the face normal
588 // will point in the wrong direction
589 lidx1 = lidx0;
590 }
591 cell_to_face_h(c1,lidx1) = f;
592
593 // Faces point from low cell index to high cell index
594 if(c0<c1){
595 face_to_cells_h(f,0) = c0;
596 face_to_localidx_h(f,0) = lidx0;
597 face_to_cells_h(f,1) = c1;
598 face_to_localidx_h(f,1) = lidx1;
599 } else {
600 face_to_cells_h(f,0) = c1;
601 face_to_localidx_h(f,0) = lidx1;
602 face_to_cells_h(f,1) = c0;
603 face_to_localidx_h(f,1) = lidx0;
604 }
605
606 // We should avoid having two virtual cells linked together.
607 TEUCHOS_ASSERT(c0<num_real_cells or c1<num_real_cells)
608
609 }
610 Kokkos::deep_copy(face_to_cells, face_to_cells_h);
611 Kokkos::deep_copy(face_to_localidx, face_to_localidx_h);
612 Kokkos::deep_copy(cell_to_face, cell_to_face_h);
613 }
614 // at this point all the data structures have been built, so now we can "do" DG.
615 // There are two of everything, an "owned" data structured corresponding to "owned"
616 // cells. And a "ghstd" data structure corresponding to ghosted cells
618 {
619 PANZER_FUNC_TIME_MONITOR_DIFF("Assign Indices",AssignIndices);
620 mesh_info.cell_to_faces = cell_to_face;
621 mesh_info.face_to_cells = face_to_cells; // faces
622 mesh_info.face_to_lidx = face_to_localidx;
623 mesh_info.subcell_dimension = space_dim;
624 mesh_info.subcell_index = -1;
625 mesh_info.has_connectivity = true;
626
627 mesh_info.num_owned_cells = owned_cells.extent(0);
628 mesh_info.num_ghstd_cells = ghost_cells.extent(0);
629 mesh_info.num_virtual_cells = virtual_cells.extent(0);
630
631 mesh_info.global_cells = PHX::View<panzer::GlobalOrdinal*>("global_cell_indices",num_total_cells);
632 mesh_info.local_cells = PHX::View<panzer::LocalOrdinal*>("local_cell_indices",num_total_cells);
633
634 Kokkos::parallel_for(num_owned_cells,KOKKOS_LAMBDA (int i) {
635 mesh_info.global_cells(i) = owned_cells(i);
636 mesh_info.local_cells(i) = i;
637 });
638
639 Kokkos::parallel_for(num_ghstd_cells,KOKKOS_LAMBDA (int i) {
640 mesh_info.global_cells(i+num_owned_cells) = ghost_cells(i);
641 mesh_info.local_cells(i+num_owned_cells) = i+num_owned_cells;
642 });
643
644 Kokkos::parallel_for(num_virtual_cells,KOKKOS_LAMBDA (int i) {
645 mesh_info.global_cells(i+num_real_cells) = virtual_cells(i);
646 mesh_info.local_cells(i+num_real_cells) = i+num_real_cells;
647 });
648
649 mesh_info.cell_nodes = PHX::View<double***>("cell_nodes",num_total_cells,nodes_per_cell,space_dim);
650
651 // Initialize coordinates to zero
652 Kokkos::deep_copy(mesh_info.cell_nodes, 0.);
653
654 Kokkos::parallel_for(num_owned_cells,KOKKOS_LAMBDA (int i) {
655 for(int j=0;j<nodes_per_cell;++j){
656 for(int k=0;k<space_dim;++k){
657 mesh_info.cell_nodes(i,j,k) = owned_nodes(i,j,k);
658 }
659 }
660 });
661
662 Kokkos::parallel_for(num_ghstd_cells,KOKKOS_LAMBDA (int i) {
663 for(int j=0;j<nodes_per_cell;++j){
664 for(int k=0;k<space_dim;++k){
665 mesh_info.cell_nodes(i+num_owned_cells,j,k) = ghstd_nodes(i,j,k);
666 }
667 }
668 });
669
670 // This will backfire at some point, but we're going to make the virtual cell have the same geometry as the cell it interfaces with
671 // This way we can define a virtual cell geometry without extruding the face outside of the domain
672 // TODO BWR Certainly, this is an issue for curved meshes
673 {
674 PANZER_FUNC_TIME_MONITOR_DIFF("Assign geometry traits",AssignGeometryTraits);
675 Kokkos::parallel_for(num_virtual_cells,KOKKOS_LAMBDA (int i) {
676 const panzer::LocalOrdinal virtual_cell = i+num_real_cells;
677 for(int local_face=0; local_face<faces_per_cell; ++local_face){
678 const panzer::LocalOrdinal face = cell_to_face(virtual_cell, local_face);
679 if(face >= 0){
680 const panzer::LocalOrdinal other_side = (face_to_cells(face, 0) == virtual_cell) ? 1 : 0;
681 const panzer::LocalOrdinal real_cell = face_to_cells(face,other_side);
682 for(int j=0;j<nodes_per_cell;++j){
683 for(int k=0;k<space_dim;++k){
684 mesh_info.cell_nodes(virtual_cell,j,k) = mesh_info.cell_nodes(real_cell,j,k);
685 }
686 }
687 break;
688 }
689 }
690 });
691
692 }
693 }
694
695 // Setup element blocks and sidesets
696 std::vector<std::string> sideset_names;
697 mesh.getSidesetNames(sideset_names);
698
699 for(const std::string & element_block_name : element_block_names){
700 PANZER_FUNC_TIME_MONITOR_DIFF("Set up setupLocalMeshBlockInfo",SetupLocalMeshBlockInfo);
701 panzer::LocalMeshBlockInfo & block_info = mesh_info.element_blocks[element_block_name];
702 setupLocalMeshBlockInfo(mesh, conn, mesh_info, element_block_name, block_info);
703 block_info.subcell_dimension = space_dim;
704 block_info.subcell_index = -1;
705 block_info.has_connectivity = true;
706
707 // Setup sidesets
708 for(const std::string & sideset_name : sideset_names){
709 PANZER_FUNC_TIME_MONITOR_DIFF("Setup LocalMeshSidesetInfo",SetupLocalMeshSidesetInfo);
710 panzer::LocalMeshSidesetInfo & sideset_info = mesh_info.sidesets[element_block_name][sideset_name];
711 setupLocalMeshSidesetInfo(mesh, conn, mesh_info, element_block_name, sideset_name, sideset_info);
712 sideset_info.subcell_dimension = space_dim;
713 sideset_info.subcell_index = -1;
714 sideset_info.has_connectivity = true;
715 }
716
717 }
718
719 return mesh_info_rcp;
720
721}
722
723}
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?
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
void getElementNodesNoResize(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 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.
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
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
PHX::View< double *** > cell_nodes
std::map< std::string, std::map< std::string, LocalMeshSidesetInfo > > sidesets
std::map< std::string, LocalMeshBlockInfo > element_blocks
Teuchos::RCP< const shards::CellTopology > cell_topology