Chaste  Release::2017.1
NodePartitioner.cpp
1 /*
2 
3 Copyright (c) 2005-2017, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 #include <cassert>
36 #include <algorithm>
37 
38 #include "Exception.hpp"
39 #include "NodePartitioner.hpp"
40 #include "PetscMatTools.hpp"
41 #include "PetscTools.hpp"
42 #include "Timer.hpp"
43 #include "TrianglesMeshReader.hpp"
44 #include "Warnings.hpp"
45 #include "petscao.h"
46 #include "petscmat.h"
47 
48 #include <parmetis.h>
49 
50 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
52  std::set<unsigned>& rNodesOwned)
53 {
54  //Note that if there is not DistributedVectorFactory in the mesh, then a dumb partition based
55  //on rMesh.GetNumNodes() is applied automatically.
56  //If there is a DistributedVectorFactory then that one will be returned
57  if (rMesh.GetDistributedVectorFactory()->GetProblemSize() != rMesh.GetNumNodes())
58  {
59  EXCEPTION("The distributed vector factory size in the mesh doesn't match the total number of nodes.");
60  }
61 
62  for (unsigned node_index = rMesh.GetDistributedVectorFactory()->GetLow();
63  node_index < rMesh.GetDistributedVectorFactory()->GetHigh();
64  node_index++)
65  {
66  rNodesOwned.insert(node_index);
67  }
68 }
69 
70 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
72  std::vector<unsigned>& rNodePermutation,
73  std::set<unsigned>& rNodesOwned,
74  std::vector<unsigned>& rProcessorsOffset)
75 {
76  assert(PetscTools::IsParallel());
77  assert(ELEMENT_DIM==2 || ELEMENT_DIM==3); // LCOV_EXCL_LINE // Metis works with triangles and tetras
78 
79  if (!PetscTools::HasParMetis()) //We must have ParMetis support compiled into Petsc
80  {
81  WARN_ONCE_ONLY("PETSc support for ParMetis is not installed."); // LCOV_EXCL_LINE
82  }
83 
84  unsigned num_nodes = rMeshReader.GetNumNodes();
85  PetscInt num_procs = PetscTools::GetNumProcs();
86  unsigned local_proc_index = PetscTools::GetMyRank();
87 
88  unsigned num_elements = rMeshReader.GetNumElements();
89  unsigned num_local_elements = num_elements / num_procs;
90  unsigned first_local_element = num_local_elements * local_proc_index;
92  {
93  // Take the excess elements
94  num_local_elements += num_elements - (num_local_elements * num_procs);
95  }
96 
98  Timer::Reset();
99 
100  /*
101  * Create PETSc matrix which will have 1 for adjacent nodes.
102  */
103  Mat connectivity_matrix;
105  PetscTools::SetupMat(connectivity_matrix, num_nodes, num_nodes, 54, PETSC_DECIDE, PETSC_DECIDE, false);
106 
107  if (!rMeshReader.IsFileFormatBinary())
108  {
109  // Advance the file pointer to the first element I own
110  for (unsigned element_index = 0; element_index < first_local_element; element_index++)
111  {
112  rMeshReader.GetNextElementData();
113  }
114  }
115 
116  // In the loop below we assume that there exist edges between any pair of nodes in an element. This is
117  // a reasonable assumption for triangles and tetrahedra. This won't be the case for quadratic simplices,
118  // squares or hexahedra (or higher order elements), leading to slightly suboptimal partitions in these
119  // cases.
120  // We allow ELEMENT_DIM smaller than SPACE_DIM in case this is a 2D mesh in
121  // a 3D space.
122  assert(SPACE_DIM >= ELEMENT_DIM);
123 
124  for (unsigned element_index = 0; element_index < num_local_elements; element_index++)
125  {
126  ElementData element_data;
127 
128  if (rMeshReader.IsFileFormatBinary())
129  {
130  element_data = rMeshReader.GetElementData(first_local_element + element_index);
131  }
132  else
133  {
134  element_data = rMeshReader.GetNextElementData();
135  }
136 
137  for (unsigned i=0; i<element_data.NodeIndices.size(); i++)
138  {
139  unsigned row = element_data.NodeIndices[i];
140  for (unsigned j=0; j<i; j++)
141  {
142  unsigned col = element_data.NodeIndices[j];
143  MatSetValue(connectivity_matrix, row, col, 1.0, INSERT_VALUES);
144  MatSetValue(connectivity_matrix, col, row, 1.0, INSERT_VALUES);
145  }
146  }
147  }
149  PetscMatTools::Finalise(connectivity_matrix);
150 
152  if (PetscTools::AmMaster())
153  {
154  Timer::PrintAndReset("Connectivity matrix assembly");
155  }
156 
157  rMeshReader.Reset();
158 
159  PetscInt connectivity_matrix_lo;
160  PetscInt connectivity_matrix_hi;
161  MatGetOwnershipRange(connectivity_matrix, &connectivity_matrix_lo, &connectivity_matrix_hi);
162 
163  unsigned num_local_nodes = connectivity_matrix_hi - connectivity_matrix_lo;
164 
165  /* PETSc MatCreateMPIAdj and parMETIS rely on adjacency arrays
166  * Named here: xadj adjncy
167  * parMETIS name: xadj adjncy (see S4.2 in parMETIS manual)
168  * PETSc name: i j (see MatCreateMPIAdj in PETSc manual)
169  *
170  * The adjacency information of all nodes is listed in the main array, adjncy. Since each node
171  * has a variable number of adjacent nodes, the array xadj is used to store the index (in adjncy) where
172  * this information starts. Since xadj[i] is the start of node i's information, xadj[i+1] marks the end.
173  *
174  *
175  */
176  MatInfo matrix_info;
177  MatGetInfo(connectivity_matrix, MAT_LOCAL, &matrix_info);
178  unsigned local_num_nz = (unsigned) matrix_info.nz_used;
179 
180  size_t size = (num_local_nodes+1)*sizeof(PetscInt);
181  void* ptr;
182  PetscMalloc(size, &ptr);
183  PetscInt* xadj = (PetscInt*) ptr;
184  size = local_num_nz*sizeof(PetscInt);
185  PetscMalloc(size, &ptr);
186  PetscInt* adjncy = (PetscInt*) ptr;
187 
188  PetscInt row_num_nz;
189  const PetscInt* column_indices;
190 
191  xadj[0]=0;
192  for (PetscInt row_global_index=connectivity_matrix_lo; row_global_index<connectivity_matrix_hi; row_global_index++)
193  {
194  MatGetRow(connectivity_matrix, row_global_index, &row_num_nz, &column_indices, PETSC_NULL);
195 
196  unsigned row_local_index = row_global_index - connectivity_matrix_lo;
197  xadj[row_local_index+1] = xadj[row_local_index] + row_num_nz;
198  for (PetscInt col_index=0; col_index<row_num_nz; col_index++)
199  {
200  adjncy[xadj[row_local_index] + col_index] = column_indices[col_index];
201  }
202 
203  MatRestoreRow(connectivity_matrix, row_global_index, &row_num_nz,&column_indices, PETSC_NULL);
204  }
205 
206  PetscTools::Destroy(connectivity_matrix);
207 
208  // Convert to an adjacency matrix
209  Mat adj_matrix;
210  MatCreateMPIAdj(PETSC_COMM_WORLD, num_local_nodes, num_nodes, xadj, adjncy, PETSC_NULL, &adj_matrix);
211 
213  if (PetscTools::AmMaster())
214  {
215  Timer::PrintAndReset("Adjacency matrix creation");
216  }
217 
218  // Get PETSc to call ParMETIS
219  MatPartitioning part;
220  MatPartitioningCreate(PETSC_COMM_WORLD, &part);
221  MatPartitioningSetAdjacency(part, adj_matrix);
222  MatPartitioningSetFromOptions(part);
223  IS new_process_numbers;
224 
225  MatPartitioningApply(part, &new_process_numbers);
226  MatPartitioningDestroy(PETSC_DESTROY_PARAM(part));
227 
229  PetscTools::Destroy(adj_matrix);
230 
232  if (PetscTools::AmMaster())
233  {
234  Timer::PrintAndReset("PETSc-ParMETIS call");
235  }
236 
237  // Figure out who owns what - processor offsets
238  PetscInt* num_nodes_per_process = new PetscInt[num_procs];
239 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
240  ISPartitioningCount(new_process_numbers, num_procs, num_nodes_per_process);
241 #else
242  ISPartitioningCount(new_process_numbers, num_nodes_per_process);
243 #endif
244 
245  rProcessorsOffset.resize(num_procs);
246  rProcessorsOffset[0] = 0;
247  for (PetscInt i=1; i<num_procs; i++)
248  {
249  rProcessorsOffset[i] = rProcessorsOffset[i-1] + num_nodes_per_process[i-1];
250  }
251  unsigned my_num_nodes = num_nodes_per_process[local_proc_index];
252  delete[] num_nodes_per_process;
253 
254  // Figure out who owns what - new node numbering
255  IS new_global_node_indices;
256  ISPartitioningToNumbering(new_process_numbers, &new_global_node_indices);
257 
258  // Index sets only give local information, we want global
259  AO ordering;
260  AOCreateBasicIS(new_global_node_indices, PETSC_NULL /* natural ordering */, &ordering);
261 
262  // The locally owned range under the new numbering
263  PetscInt* local_range = new PetscInt[my_num_nodes];
264  for (unsigned i=0; i<my_num_nodes; i++)
265  {
266  local_range[i] = rProcessorsOffset[local_proc_index] + i;
267  }
268  AOApplicationToPetsc(ordering, my_num_nodes, local_range);
269  //AOView(ordering, PETSC_VIEWER_STDOUT_WORLD);
270 
271  // Fill in rNodesOwned
272  rNodesOwned.insert(local_range, local_range + my_num_nodes);
273  delete[] local_range;
274 
275  // Once we know the offsets we can compute the permutation vector
276  PetscInt* global_range = new PetscInt[num_nodes];
277  for (unsigned i=0; i<num_nodes; i++)
278  {
279  global_range[i] = i;
280  }
281  AOPetscToApplication(ordering, num_nodes, global_range);
282 
283  rNodePermutation.resize(num_nodes);
284  std::copy(global_range, global_range+num_nodes, rNodePermutation.begin());
285  delete[] global_range;
286 
287  AODestroy(PETSC_DESTROY_PARAM(ordering));
288  ISDestroy(PETSC_DESTROY_PARAM(new_process_numbers));
289  ISDestroy(PETSC_DESTROY_PARAM(new_global_node_indices));
290 
292  if (PetscTools::AmMaster())
293  {
294  Timer::PrintAndReset("PETSc-ParMETIS output manipulation");
295  }
296 }
297 
298 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
300  std::vector<unsigned>& rNodePermutation,
301  std::set<unsigned>& rNodesOwned,
302  std::vector<unsigned>& rProcessorsOffset,
303  ChasteCuboid<SPACE_DIM>* pRegion)
304 {
306 
307  unsigned num_nodes = rMeshReader.GetNumNodes();
308 
309  // Work out where each node should lie.
310  std::vector<unsigned> node_ownership(num_nodes, 0);
311  std::vector<unsigned> node_index_ownership(num_nodes, UNSIGNED_UNSET);
312 
313  for (unsigned node=0; node < num_nodes; node++)
314  {
315  std::vector<double> location = rMeshReader.GetNextNode();
316 
317  // Make sure it is the correct size
318  assert(location.size() == SPACE_DIM);
319 
320  // Establish whether it lies in the domain. ChasteCuboid::DoesContain is
321  // insufficient for this as it treats all boundaries as open.
322  const ChastePoint<SPACE_DIM>& lower = pRegion->rGetLowerCorner();
323  const ChastePoint<SPACE_DIM>& upper = pRegion->rGetUpperCorner();
324 
325  bool does_contain = true;
326 
327  for (unsigned d=0; d<SPACE_DIM; d++)
328  {
329  bool boundary_check;
330  boundary_check = ((location[d] > lower[d]) || sqrt((location[d]-lower[d])*(location[d]-lower[d])) < DBL_EPSILON);
331  does_contain = (does_contain && boundary_check && (location[d] < upper[d]));
332  }
333 
334  if (does_contain)
335  {
336  node_ownership[node] = 1;
337  rNodesOwned.insert(node);
338  node_index_ownership[node] = PetscTools::GetMyRank();
339  }
340  }
341 
342  // Make sure each node will be owned by exactly on process
343  std::vector<unsigned> global_ownership(num_nodes, 0);
344  std::vector<unsigned> global_index_ownership(num_nodes, UNSIGNED_UNSET);
345 
346  MPI_Allreduce(&node_ownership[0], &global_ownership[0], num_nodes, MPI_UNSIGNED, MPI_LXOR, PETSC_COMM_WORLD);
347  for (unsigned i=0; i<num_nodes; i++)
348  {
349  if (global_ownership[i] == 0)
350  {
351  EXCEPTION("A node is either not in geometric region, or the regions are not disjoint.");
352  }
353  }
354 
355  // Create the permutation and offset vectors.
356  MPI_Allreduce(&node_index_ownership[0], &global_index_ownership[0], num_nodes, MPI_UNSIGNED, MPI_MIN, PETSC_COMM_WORLD);
357 
358  for (unsigned proc=0; proc<PetscTools::GetNumProcs(); proc++)
359  {
360  rProcessorsOffset.push_back(rNodePermutation.size());
361  for (unsigned node=0; node<num_nodes; node++)
362  {
363  if (global_index_ownership[node] == proc)
364  {
365  rNodePermutation.push_back(node);
366  }
367  }
368  }
369  assert(rNodePermutation.size() == num_nodes);
370 }
371 
372 // Explicit instantiation
373 template class NodePartitioner<1,1>;
374 template class NodePartitioner<1,2>;
375 template class NodePartitioner<1,3>;
376 template class NodePartitioner<2,2>;
377 template class NodePartitioner<2,3>;
378 template class NodePartitioner<3,3>;
virtual DistributedVectorFactory * GetDistributedVectorFactory()
virtual ElementData GetNextElementData()=0
virtual ElementData GetElementData(unsigned index)
static void Barrier(const std::string callerId="")
Definition: PetscTools.cpp:134
#define PETSC_DESTROY_PARAM(x)
Definition: PetscTools.hpp:70
#define EXCEPTION(message)
Definition: Exception.hpp:143
static bool AmMaster()
Definition: PetscTools.cpp:120
static void DumbPartitioning(AbstractMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, std::set< unsigned > &rNodesOwned)
static void GeometricPartitioning(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader, std::vector< unsigned > &rNodePermutation, std::set< unsigned > &rNodesOwned, std::vector< unsigned > &rProcessorsOffset, ChasteCuboid< SPACE_DIM > *pRegion)
virtual unsigned GetNumNodes() const
static void PrintAndReset(std::string message)
Definition: Timer.cpp:70
std::vector< unsigned > NodeIndices
const ChastePoint< SPACE_DIM > & rGetUpperCorner() const
static void PetscMatrixPartitioning(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader, std::vector< unsigned > &rNodePermutation, std::set< unsigned > &rNodesOwned, std::vector< unsigned > &rProcessorsOffset)
virtual unsigned GetNumElements() const =0
virtual void Reset()=0
static bool HasParMetis()
Definition: PetscTools.cpp:445
const unsigned UNSIGNED_UNSET
Definition: Exception.hpp:52
static void SetupMat(Mat &rMat, int numRows, int numColumns, unsigned rowPreallocation, int numLocalRows=PETSC_DECIDE, int numLocalColumns=PETSC_DECIDE, bool ignoreOffProcEntries=true, bool newAllocationError=true)
Definition: PetscTools.cpp:268
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:352
static bool IsParallel()
Definition: PetscTools.cpp:97
virtual bool IsFileFormatBinary()
static void Finalise(Mat matrix)
static bool AmTopMost()
Definition: PetscTools.cpp:126
static void Reset()
Definition: Timer.cpp:44
static unsigned GetNumProcs()
Definition: PetscTools.cpp:108
virtual std::vector< double > GetNextNode()=0
static unsigned GetMyRank()
Definition: PetscTools.cpp:114
virtual unsigned GetNumNodes() const =0
const ChastePoint< SPACE_DIM > & rGetLowerCorner() const