Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_LocalPartitioningUtilities.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
44 
45 #include "Teuchos_RCP.hpp"
46 #include "Teuchos_Comm.hpp"
47 #include "Teuchos_Assert.hpp"
48 
49 #include "Panzer_Workset_Builder.hpp"
51 
52 #include "Phalanx_KokkosDeviceTypes.hpp"
53 
54 #include <unordered_set>
55 #include <unordered_map>
56 
57 namespace panzer
58 {
59 
60 namespace partitioning_utilities
61 {
62 
63 void
65  const std::vector<panzer::LocalOrdinal> & owned_parent_cells,
66  panzer::LocalMeshInfoBase & sub_info)
67 {
68  using GO = panzer::GlobalOrdinal;
69  using LO = panzer::LocalOrdinal;
70 
71  PANZER_FUNC_TIME_MONITOR_DIFF("panzer::partitioning_utilities::setupSubLocalMeshInfo",setupSLMI);
72  // The goal of this function is to fill a LocalMeshInfoBase (sub_info) with
73  // a subset of cells from a given parent LocalMeshInfoBase (parent_info)
74 
75  // Note: owned_parent_cells are the owned cells for sub_info in the parent_info's indexing scheme
76  // We need to generate sub_info's ghosts and figure out the virtual cells
77 
78  // Note: We will only handle a single ghost layer
79 
80  // Note: We assume owned_parent_cells are owned cells of the parent
81  // i.e. owned_parent_indexes cannot refer to ghost or virtual cells in parent_info
82 
83  // Note: This function works with inter-face connectivity. NOT node connectivity.
84 
85  const int num_owned_cells = owned_parent_cells.size();
86  TEUCHOS_TEST_FOR_EXCEPT_MSG(num_owned_cells == 0, "panzer::partitioning_utilities::setupSubLocalMeshInfo : Input parent subcells must exist (owned_parent_cells)");
87 
88  const int num_parent_owned_cells = parent_info.num_owned_cells;
89  TEUCHOS_TEST_FOR_EXCEPT_MSG(num_parent_owned_cells <= 0, "panzer::partitioning_utilities::setupSubLocalMeshInfo : Input parent info must contain owned cells");
90 
91  const int num_parent_ghstd_cells = parent_info.num_ghstd_cells;
92 
93  const int num_parent_total_cells = parent_info.num_owned_cells + parent_info.num_ghstd_cells + parent_info.num_virtual_cells;
94 
95  // Just as a precaution, make sure the parent_info is setup properly
96  TEUCHOS_ASSERT(static_cast<int>(parent_info.cell_to_faces.extent(0)) == num_parent_total_cells);
97  const int num_faces_per_cell = parent_info.cell_to_faces.extent(1);
98 
99  // The first thing to do is construct a vector containing the parent cell indexes of all
100  // owned, ghstd, and virtual cells
101  std::vector<LO> ghstd_parent_cells;
102  std::vector<LO> virtual_parent_cells;
103  {
104  PANZER_FUNC_TIME_MONITOR_DIFF("Construct parent cell vector",ParentCell);
105  // We grab all of the owned cells and put their global indexes into sub_info
106  // We also put all of the owned cell indexes in the parent's indexing scheme into a set to use for lookups
107  std::unordered_set<LO> owned_parent_cells_set(owned_parent_cells.begin(), owned_parent_cells.end());
108 
109  // We need to create a list of ghstd and virtual cells
110  // We do this by running through sub_cell_indexes
111  // and looking at the neighbors to find neighbors that are not owned
112 
113  // Virtual cells are defined as cells with indexes outside of the range of owned_cells and ghstd_cells
114  const int virtual_parent_cell_offset = num_parent_owned_cells + num_parent_ghstd_cells;
115 
116  std::unordered_set<LO> ghstd_parent_cells_set;
117  std::unordered_set<LO> virtual_parent_cells_set;
118  for(int i=0;i<num_owned_cells;++i){
119  const LO parent_cell_index = owned_parent_cells[i];
120  for(int local_face_index=0;local_face_index<num_faces_per_cell;++local_face_index){
121  const LO parent_face = parent_info.cell_to_faces(parent_cell_index, local_face_index);
122 
123  // Sidesets can have owned cells that border the edge of the domain (i.e. parent_face == -1)
124  // If we are at the edge of the domain, we can ignore this face.
125  if(parent_face < 0)
126  continue;
127 
128  // Find the side index for neighbor cell with respect to the face
129  const LO neighbor_parent_side = (parent_info.face_to_cells(parent_face,0) == parent_cell_index) ? 1 : 0;
130 
131  // Get the neighbor cell index in the parent's indexing scheme
132  const LO neighbor_parent_cell = parent_info.face_to_cells(parent_face, neighbor_parent_side);
133 
134  // If the face exists, then the neighbor should exist
135  TEUCHOS_ASSERT(neighbor_parent_cell >= 0);
136 
137  // We can easily check if this is a virtual cell
138  if(neighbor_parent_cell >= virtual_parent_cell_offset){
139  virtual_parent_cells_set.insert(neighbor_parent_cell);
140  } else if(neighbor_parent_cell >= num_parent_owned_cells){
141  // This is a quick check for a ghost cell
142  // This branch only exists to cut down on the number of times the next branch (much slower) is called
143  ghstd_parent_cells_set.insert(neighbor_parent_cell);
144  } else {
145  // There is still potential for this to be a ghost cell with respect to 'our' cells
146  // The only way to check this is with a super slow lookup call
147  if(owned_parent_cells_set.find(neighbor_parent_cell) == owned_parent_cells_set.end()){
148  // The neighbor cell is not owned by 'us', therefore it is a ghost
149  ghstd_parent_cells_set.insert(neighbor_parent_cell);
150  }
151  }
152  }
153  }
154 
155  // We now have a list of the owned, ghstd, and virtual cells in the parent's indexing scheme.
156  // We will take the 'unordered_set's ordering for the the sub-indexing scheme
157 
158  ghstd_parent_cells.assign(ghstd_parent_cells_set.begin(), ghstd_parent_cells_set.end());
159  virtual_parent_cells.assign(virtual_parent_cells_set.begin(), virtual_parent_cells_set.end());
160 
161  }
162 
163  const int num_ghstd_cells = ghstd_parent_cells.size();
164  const int num_virtual_cells = virtual_parent_cells.size();
165  const int num_real_cells = num_owned_cells + num_ghstd_cells;
166  const int num_total_cells = num_real_cells + num_virtual_cells;
167 
168  std::vector<std::pair<LO, LO> > all_parent_cells(num_total_cells);
169  for (std::size_t i=0; i< owned_parent_cells.size(); ++i)
170  all_parent_cells[i] = std::pair<LO, LO>(owned_parent_cells[i], i);
171 
172  for (std::size_t i=0; i< ghstd_parent_cells.size(); ++i) {
173  LO insert = owned_parent_cells.size()+i;
174  all_parent_cells[insert] = std::pair<LO, LO>(ghstd_parent_cells[i], insert);
175  }
176 
177  for (std::size_t i=0; i< virtual_parent_cells.size(); ++i) {
178  LO insert = owned_parent_cells.size()+ ghstd_parent_cells.size()+i;
179  all_parent_cells[insert] = std::pair<LO, LO>(virtual_parent_cells[i], insert);
180  }
181 
182  sub_info.num_owned_cells = owned_parent_cells.size();
183  sub_info.num_ghstd_cells = ghstd_parent_cells.size();
184  sub_info.num_virtual_cells = virtual_parent_cells.size();
185 
186  // We now have the indexing order for our sub_info
187 
188  // Just as a precaution, make sure the parent_info is setup properly
189  TEUCHOS_ASSERT(static_cast<int>(parent_info.cell_vertices.extent(0)) == num_parent_total_cells);
190  TEUCHOS_ASSERT(static_cast<int>(parent_info.local_cells.extent(0)) == num_parent_total_cells);
191  TEUCHOS_ASSERT(static_cast<int>(parent_info.global_cells.extent(0)) == num_parent_total_cells);
192 
193  const int num_vertices_per_cell = parent_info.cell_vertices.extent(1);
194  const int num_dims = parent_info.cell_vertices.extent(2);
195 
196  // Fill owned, ghstd, and virtual cells: global indexes, local indexes and vertices
197  sub_info.global_cells = Kokkos::View<GO*>("global_cells", num_total_cells);
198  sub_info.local_cells = Kokkos::View<LO*>("local_cells", num_total_cells);
199  sub_info.cell_vertices = Kokkos::View<double***, PHX::Device>("cell_vertices", num_total_cells, num_vertices_per_cell, num_dims);
200  for(int cell=0;cell<num_total_cells;++cell){
201  const LO parent_cell = all_parent_cells[cell].first;
202  sub_info.global_cells(cell) = parent_info.global_cells(parent_cell);
203  sub_info.local_cells(cell) = parent_info.local_cells(parent_cell);
204  for(int vertex=0;vertex<num_vertices_per_cell;++vertex){
205  for(int dim=0;dim<num_dims;++dim){
206  sub_info.cell_vertices(cell,vertex,dim) = parent_info.cell_vertices(parent_cell,vertex,dim);
207  }
208  }
209  }
210 
211  // Now for the difficult part
212 
213  // We need to create a new face indexing scheme from the old face indexing scheme
214 
215  // Create an auxiliary list with all cells - note this preserves indexing
216 
217  struct face_t{
218  face_t(LO c0, LO c1, LO sc0, LO sc1)
219  {
220  cell_0=c0;
221  cell_1=c1;
222  subcell_index_0=sc0;
223  subcell_index_1=sc1;
224  }
225  LO cell_0;
226  LO cell_1;
227  LO subcell_index_0;
228  LO subcell_index_1;
229  };
230 
231 
232  // First create the faces
233  std::vector<face_t> faces;
234  {
235  PANZER_FUNC_TIME_MONITOR_DIFF("Create faces",CreateFaces);
236  // faces_set: cell_0, subcell_index_0, cell_1, subcell_index_1
237  std::unordered_map<LO,std::unordered_map<LO, std::pair<LO,LO> > > faces_set;
238 
239  std::sort(all_parent_cells.begin(), all_parent_cells.end());
240 
241  for(int owned_cell=0;owned_cell<num_owned_cells;++owned_cell){
242  const LO owned_parent_cell = owned_parent_cells[owned_cell];
243  for(int local_face=0;local_face<num_faces_per_cell;++local_face){
244  const LO parent_face = parent_info.cell_to_faces(owned_parent_cell,local_face);
245 
246  // Skip faces at the edge of the domain
247  if(parent_face<0)
248  continue;
249 
250  // Get the cell on the other side of the face
251  const LO neighbor_side = (parent_info.face_to_cells(parent_face,0) == owned_parent_cell) ? 1 : 0;
252 
253  const LO neighbor_parent_cell = parent_info.face_to_cells(parent_face, neighbor_side);
254  const LO neighbor_subcell_index = parent_info.face_to_lidx(parent_face, neighbor_side);
255 
256  // Convert parent cell index into sub cell index
257  std::pair<LO, LO> search_point(neighbor_parent_cell, 0);
258  auto itr = std::lower_bound(all_parent_cells.begin(), all_parent_cells.end(), search_point);
259 
260  TEUCHOS_TEST_FOR_EXCEPT_MSG(itr == all_parent_cells.end(), "panzer_stk::setupSubLocalMeshInfo : Neighbor cell was not found in owned, ghosted, or virtual cells");
261 
262  const LO neighbor_cell = itr->second;
263 
264  LO cell_0, cell_1, subcell_index_0, subcell_index_1;
265  if(owned_cell < neighbor_cell){
266  cell_0 = owned_cell;
267  subcell_index_0 = local_face;
268  cell_1 = neighbor_cell;
269  subcell_index_1 = neighbor_subcell_index;
270  } else {
271  cell_1 = owned_cell;
272  subcell_index_1 = local_face;
273  cell_0 = neighbor_cell;
274  subcell_index_0 = neighbor_subcell_index;
275  }
276 
277  // Add this interface to the set of faces - smaller cell index is 'left' (or '0') side of face
278  faces_set[cell_0][subcell_index_0] = std::pair<LO,LO>(cell_1, subcell_index_1);
279  }
280  }
281 
282  for(const auto & cell_pair : faces_set){
283  const LO cell_0 = cell_pair.first;
284  for(const auto & subcell_pair : cell_pair.second){
285  const LO subcell_index_0 = subcell_pair.first;
286  const LO cell_1 = subcell_pair.second.first;
287  const LO subcell_index_1 = subcell_pair.second.second;
288  faces.push_back(face_t(cell_0,cell_1,subcell_index_0,subcell_index_1));
289  }
290  }
291  }
292 
293  const int num_faces = faces.size();
294 
295  sub_info.face_to_cells = Kokkos::View<LO*[2]>("face_to_cells", num_faces);
296  sub_info.face_to_lidx = Kokkos::View<LO*[2]>("face_to_lidx", num_faces);
297  sub_info.cell_to_faces = Kokkos::View<LO**>("cell_to_faces", num_total_cells, num_faces_per_cell);
298 
299  // Default the system with invalid cell index
300  Kokkos::deep_copy(sub_info.cell_to_faces, -1);
301 
302  for(int face_index=0;face_index<num_faces;++face_index){
303  const face_t & face = faces[face_index];
304 
305  sub_info.face_to_cells(face_index,0) = face.cell_0;
306  sub_info.face_to_cells(face_index,1) = face.cell_1;
307 
308  sub_info.cell_to_faces(face.cell_0,face.subcell_index_0) = face_index;
309  sub_info.cell_to_faces(face.cell_1,face.subcell_index_1) = face_index;
310 
311  sub_info.face_to_lidx(face_index,0) = face.subcell_index_0;
312  sub_info.face_to_lidx(face_index,1) = face.subcell_index_1;
313 
314  }
315 
316  // Complete.
317 
318 }
319 
320 void
322  const int splitting_size,
323  std::vector<panzer::LocalMeshPartition> & partitions)
324 {
325 
326  using LO = panzer::LocalOrdinal;
327 
328  // Make sure the splitting size makes sense
329  TEUCHOS_ASSERT((splitting_size > 0) or (splitting_size == WorksetSizeType::ALL_ELEMENTS));
330 
331  // Default partition size
332  const LO base_partition_size = std::min(mesh_info.num_owned_cells, (splitting_size > 0) ? splitting_size : mesh_info.num_owned_cells);
333 
334  // Cells to partition
335  std::vector<LO> partition_cells;
336  partition_cells.resize(base_partition_size);
337 
338  // Create the partitions
339  LO cell_count = 0;
340  while(cell_count < mesh_info.num_owned_cells){
341 
342  LO partition_size = base_partition_size;
343  if(cell_count + partition_size > mesh_info.num_owned_cells)
344  partition_size = mesh_info.num_owned_cells - cell_count;
345 
346  // Error check for a null partition - this should never happen by design
347  TEUCHOS_ASSERT(partition_size != 0);
348 
349  // In the final partition, we need to reduce the size of partition_cells
350  if(partition_size != base_partition_size)
351  partition_cells.resize(partition_size);
352 
353  // Set the partition indexes - not really a partition, just a chunk of cells
354  for(LO i=0; i<partition_size; ++i)
355  partition_cells[i] = cell_count+i;
356 
357  // Create an empty partition
358  partitions.push_back(panzer::LocalMeshPartition());
359 
360  // Fill the empty partition
361  partitioning_utilities::setupSubLocalMeshInfo(mesh_info,partition_cells,partitions.back());
362 
363  // Update the cell count
364  cell_count += partition_size;
365  }
366 
367 }
368 
369 }
370 
371 void
373  const panzer::WorksetDescriptor & description,
374  std::vector<panzer::LocalMeshPartition> & partitions)
375 {
376  // We have to make sure that the partitioning is possible
378  TEUCHOS_ASSERT(description.getWorksetSize() != 0);
379 
380  // This could just return, but it would be difficult to debug why no partitions were returned
381  TEUCHOS_ASSERT(description.requiresPartitioning());
382 
383  const std::string & element_block_name = description.getElementBlock();
384 
385  // We have two processes for in case this is a sideset or element block
386  if(description.useSideset()){
387 
388  // If the element block doesn't exist, there are no partitions to create
389  if(mesh_info.sidesets.find(element_block_name) == mesh_info.sidesets.end())
390  return;
391  const auto & sideset_map = mesh_info.sidesets.at(element_block_name);
392 
393  const std::string & sideset_name = description.getSideset();
394 
395  // If the sideset doesn't exist, there are no partitions to create
396  if(sideset_map.find(sideset_name) == sideset_map.end())
397  return;
398 
399  const panzer::LocalMeshSidesetInfo & sideset_info = sideset_map.at(sideset_name);
400 
401  // Partitioning is not important for sidesets
402  panzer::partitioning_utilities::splitMeshInfo(sideset_info, description.getWorksetSize(), partitions);
403 
404  for(auto & partition : partitions){
405  partition.sideset_name = sideset_name;
406  partition.element_block_name = element_block_name;
407  partition.cell_topology = sideset_info.cell_topology;
408  }
409 
410  } else {
411 
412  // If the element block doesn't exist, there are no partitions to create
413  if(mesh_info.element_blocks.find(element_block_name) == mesh_info.element_blocks.end())
414  return;
415 
416  // Grab the element block we're interested in
417  const panzer::LocalMeshBlockInfo & block_info = mesh_info.element_blocks.at(element_block_name);
418 
420  // We only have one partition describing the entire local mesh
421  panzer::partitioning_utilities::splitMeshInfo(block_info, -1, partitions);
422  } else {
423  // We need to partition local mesh
424 
425  // FIXME: Until the above function is fixed, we will use this hack - this will lead to horrible partitions
426  panzer::partitioning_utilities::splitMeshInfo(block_info, description.getWorksetSize(), partitions);
427 
428  }
429 
430  for(auto & partition : partitions){
431  partition.element_block_name = element_block_name;
432  partition.cell_topology = block_info.cell_topology;
433  }
434  }
435 
436 }
437 
438 }
Backwards compatibility mode that ignores the worksetSize in the WorksetDescriptor.
std::map< std::string, std::map< std::string, LocalMeshSidesetInfo > > sidesets
Kokkos::View< double ***, PHX::Device > cell_vertices
const std::string & getSideset(const int block=0) const
Get sideset name.
Kokkos::View< panzer::LocalOrdinal * > local_cells
Teuchos::RCP< const shards::CellTopology > cell_topology
const std::string & getElementBlock(const int block=0) const
Get element block name.
void generateLocalMeshPartitions(const panzer::LocalMeshInfo &mesh_info, const panzer::WorksetDescriptor &description, std::vector< panzer::LocalMeshPartition > &partitions)
bool useSideset() const
This descriptor is for a side set.
Kokkos::View< panzer::GlobalOrdinal * > global_cells
panzer::LocalOrdinal num_owned_cells
Kokkos::View< panzer::LocalOrdinal *[2]> face_to_cells
panzer::LocalOrdinal num_virtual_cells
Teuchos::RCP< const shards::CellTopology > cell_topology
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
panzer::LocalOrdinal num_ghstd_cells
std::map< std::string, LocalMeshBlockInfo > element_blocks
Kokkos::View< panzer::LocalOrdinal *[2]> face_to_lidx
Kokkos::View< panzer::LocalOrdinal ** > cell_to_faces
void splitMeshInfo(const panzer::LocalMeshInfoBase &mesh_info, const int splitting_size, std::vector< panzer::LocalMeshPartition > &partitions)
#define TEUCHOS_ASSERT(assertion_test)
int getWorksetSize() const
Get the requested workset size (default -2 (workset size is set elsewhere), -1 (largest possible work...
void setupSubLocalMeshInfo(const panzer::LocalMeshInfoBase &parent_info, const std::vector< panzer::LocalOrdinal > &owned_parent_cells, panzer::LocalMeshInfoBase &sub_info)
bool requiresPartitioning() const
Do we need to partition the local mesh prior to generating worksets.
Workset size is set to the total number of local elements in the MPI process.