Chaste  Release::3.4
DistributedBoxCollection.cpp
1 /*
2 
3 Copyright (c) 2005-2016, 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 "DistributedBoxCollection.hpp"
36 #include "Exception.hpp"
37 #include "MathsCustomFunctions.hpp"
38 #include "Warnings.hpp"
39 
40 // Static member for "fudge factor" is instantiated here
41 template<unsigned DIM>
42 const double DistributedBoxCollection<DIM>::msFudge = 5e-14;
43 
44 template<unsigned DIM>
45 DistributedBoxCollection<DIM>::DistributedBoxCollection(double boxWidth, c_vector<double, 2*DIM> domainSize, bool isPeriodicInX, int localRows)
46  : mBoxWidth(boxWidth),
47  mIsPeriodicInX(isPeriodicInX),
48  mAreLocalBoxesSet(false),
49  mCalculateNodeNeighbours(true)
50 {
51  // Periodicity only works in 2d and in serial.
52  if (isPeriodicInX)
53  {
54  assert(DIM==2 && PetscTools::IsSequential());
55  }
56 
57  // If the domain size is not 'divisible' (i.e. fmod(width, box_size) > 0.0) we swell the domain to enforce this.
58  for (unsigned i=0; i<DIM; i++)
59  {
60  double r = fmod((domainSize[2*i+1]-domainSize[2*i]), boxWidth);
61  if (r > 0.0)
62  {
63  domainSize[2*i+1] += boxWidth - r;
64  }
65  }
66 
67  mDomainSize = domainSize;
68 
69  // Calculate the number of boxes in each direction.
70  mNumBoxesEachDirection = scalar_vector<unsigned>(DIM, 0u);
71 
72  for (unsigned i=0; i<DIM; i++)
73  {
74  double counter = mDomainSize(2*i);
75  while (counter + msFudge < mDomainSize(2*i+1))
76  {
78  counter += mBoxWidth;
79  }
80  }
81 
82  // Make sure there are enough boxes for the number of processes.
84  {
85  WARNING("There are more processes than convenient for the domain/mesh/box size. The domain size has been swollen.")
88  }
89 
90  // Make a distributed vectory factory to split the rows of boxes between processes.
92 
93  // Calculate how many boxes in a row / face. A useful piece of data in the class.
94  mNumBoxesInAFace = 1;
95  for (unsigned i=1; i<DIM; i++)
96  {
98  }
99 
101 
104 
105  /*
106  * The location of the Boxes doesn't matter as it isn't actually used so we don't bother to work it out.
107  * The reason it isn't used is because this class works out which box a node lies in without refernece to actual
108  * box locations, only their global index within the collection.
109  */
110  c_vector<double, 2*DIM> arbitrary_location;
111  for (unsigned i=0; i<num_boxes; i++)
112  {
113  Box<DIM> new_box(arbitrary_location);
114  mBoxes.push_back(new_box);
115 
116  unsigned global_index = mMinBoxIndex + i;
117  mBoxesMapping[global_index] = mBoxes.size() - 1;
118  }
119 
121 }
122 
123 template<unsigned DIM>
125 {
126  delete mpDistributedBoxStackFactory;
127 }
128 
129 template<unsigned DIM>
131 {
132  for (unsigned i=0; i<mBoxes.size(); i++)
133  {
134  mBoxes[i].ClearNodes();
135  }
136  for (unsigned i=0; i<mHaloBoxes.size(); i++)
137  {
138  mHaloBoxes[i].ClearNodes();
139  }
140 }
141 
142 template<unsigned DIM>
144 {
145  // Get top-most and bottom-most value of Distributed Box Stack.
146  unsigned Hi = mpDistributedBoxStackFactory->GetHigh();
147  unsigned Lo = mpDistributedBoxStackFactory->GetLow();
148 
149  // If I am not the top-most process, add halo structures to the right.
150  if (!PetscTools::AmTopMost())
151  {
152  for (unsigned i=0; i < mNumBoxesInAFace; i++)
153  {
154  c_vector<double, 2*DIM> arbitrary_location; // See comment in constructor.
155  Box<DIM> new_box(arbitrary_location);
156  mHaloBoxes.push_back(new_box);
157 
158  unsigned global_index = Hi * mNumBoxesInAFace + i;
159  mHaloBoxesMapping[global_index] = mHaloBoxes.size()-1;
160  mHalosRight.push_back(global_index - mNumBoxesInAFace);
161  }
162  }
163 
164  // If I am not the bottom-most process, add halo structures to the left.
165  if (!PetscTools::AmMaster())
166  {
167  for (unsigned i=0; i< mNumBoxesInAFace; i++)
168  {
169  c_vector<double, 2*DIM> arbitrary_location; // See comment in constructor.
170  Box<DIM> new_box(arbitrary_location);
171  mHaloBoxes.push_back(new_box);
172 
173  unsigned global_index = (Lo - 1) * mNumBoxesInAFace + i;
174  mHaloBoxesMapping[global_index] = mHaloBoxes.size() - 1;
175 
176  mHalosLeft.push_back(global_index + mNumBoxesInAFace);
177  }
178  }
179 }
180 
181 template<unsigned DIM>
183 {
184  mHaloNodesLeft.clear();
185  for (unsigned i=0; i<mHalosLeft.size(); i++)
186  {
187  for (typename std::set<Node<DIM>* >::iterator iter=this->rGetBox(mHalosLeft[i]).rGetNodesContained().begin();
188  iter!=this->rGetBox(mHalosLeft[i]).rGetNodesContained().end();
189  iter++)
190  {
191  mHaloNodesLeft.push_back((*iter)->GetIndex());
192  }
193  }
194 
195  // Send right
196  mHaloNodesRight.clear();
197  for (unsigned i=0; i<mHalosRight.size(); i++)
198  {
199  for (typename std::set<Node<DIM>* >::iterator iter=this->rGetBox(mHalosRight[i]).rGetNodesContained().begin();
200  iter!=this->rGetBox(mHalosRight[i]).rGetNodesContained().end();
201  iter++)
202  {
203  mHaloNodesRight.push_back((*iter)->GetIndex());
204  }
205  }
206 }
207 
208 template<unsigned DIM>
210 {
211  return mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow();
212 }
213 
214 template<unsigned DIM>
216 {
217  return (!(globalIndex<mMinBoxIndex) && !(mMaxBoxIndex<globalIndex));
218 }
219 
220 template<unsigned DIM>
222 {
223  bool is_halo_right = ((globalIndex > mMaxBoxIndex) && !(globalIndex > mMaxBoxIndex + mNumBoxesInAFace));
224  bool is_halo_left = ((globalIndex < mMinBoxIndex) && !(globalIndex < mMinBoxIndex - mNumBoxesInAFace));
225 
226  return (PetscTools::IsParallel() && (is_halo_right || is_halo_left));
227 }
228 
229 template<unsigned DIM>
231 {
232  bool is_on_boundary = !(globalIndex < mMaxBoxIndex - mNumBoxesInAFace) || (globalIndex < mMinBoxIndex + mNumBoxesInAFace);
233 
234  return (PetscTools::IsSequential() || !(is_on_boundary));
235 }
236 
237 template<unsigned DIM>
238 unsigned DistributedBoxCollection<DIM>::CalculateGlobalIndex(c_vector<unsigned, DIM> coordinateIndices)
239 {
240  unsigned containing_box_index = 0;
241  for (unsigned i=0; i<DIM; i++)
242  {
243  unsigned temp = 1;
244  for (unsigned j=0; j<i; j++)
245  {
246  temp *= mNumBoxesEachDirection(j);
247  }
248  containing_box_index += temp*coordinateIndices[i];
249  }
250 
251  return containing_box_index;
252 }
253 
254 template<unsigned DIM>
256 {
257  // Get the location of the node
258  c_vector<double, DIM> location = pNode->rGetLocation();
259  return CalculateContainingBox(location);
260 }
261 
262 
263 template<unsigned DIM>
264 unsigned DistributedBoxCollection<DIM>::CalculateContainingBox(c_vector<double, DIM>& rLocation)
265 {
266  // The node must lie inside the boundary of the box collection
267  for (unsigned i=0; i<DIM; i++)
268  {
269  if ( (rLocation[i] < mDomainSize(2*i)) || !(rLocation[i] < mDomainSize(2*i+1)) )
270  {
271  EXCEPTION("The point provided is outside all of the boxes");
272  }
273  }
274 
275  // Compute the containing box index in each dimension
276  c_vector<unsigned, DIM> containing_box_indices = scalar_vector<unsigned>(DIM, 0u);
277  for (unsigned i=0; i<DIM; i++)
278  {
279  double box_counter = mDomainSize(2*i);
280  while (!((box_counter + mBoxWidth) > rLocation[i] + msFudge))
281  {
282  containing_box_indices[i]++;
283  box_counter += mBoxWidth;
284  }
285  }
286 
287  // Use these to compute the index of the containing box
288  unsigned containing_box_index = 0;
289  for (unsigned i=0; i<DIM; i++)
290  {
291  unsigned temp = 1;
292  for (unsigned j=0; j<i; j++)
293  {
294  temp *= mNumBoxesEachDirection(j);
295  }
296  containing_box_index += temp*containing_box_indices[i];
297  }
298 
299  // This index must be less than the total number of boxes
300  assert(containing_box_index < mNumBoxes);
301 
302  return containing_box_index;
303 }
304 
305 template<unsigned DIM>
306 c_vector<unsigned, DIM> DistributedBoxCollection<DIM>::CalculateCoordinateIndices(unsigned globalIndex)
307 {
308  c_vector<unsigned, DIM> indices;
309 
310  switch(DIM)
311  {
312  case 1:
313  {
314  indices[0]=globalIndex;
315  break;
316  }
317  case 2:
318  {
319  unsigned remainder=globalIndex % mNumBoxesEachDirection(0);
320  indices[0]=remainder;
321  indices[1]=(unsigned)(globalIndex/mNumBoxesEachDirection(0));
322  break;
323  }
324 
325  case 3:
326  {
327  unsigned remainder1=globalIndex % (mNumBoxesEachDirection(0)*mNumBoxesEachDirection(1));
328  unsigned remainder2=remainder1 % mNumBoxesEachDirection(0);
329  indices[0]=remainder2;
330  indices[1]=((globalIndex-indices[0])/mNumBoxesEachDirection(0))%mNumBoxesEachDirection(1);
331  indices[2]=((globalIndex-indices[0]-mNumBoxesEachDirection(0)*indices[1])/(mNumBoxesEachDirection(0)*mNumBoxesEachDirection(1)));
332  break;
333  }
334  default:
335  {
337  }
338  }
339 
340  return indices;
341 }
342 
343 template<unsigned DIM>
345 {
346  assert(!(boxIndex<mMinBoxIndex) && !(mMaxBoxIndex<boxIndex));
347  return mBoxes[boxIndex-mMinBoxIndex];
348 }
349 
350 template<unsigned DIM>
352 {
353  assert(GetHaloBoxOwnership(boxIndex));
354 
355  unsigned local_index = mHaloBoxesMapping.find(boxIndex)->second;
356 
357  return mHaloBoxes[local_index];
358 }
359 
360 template<unsigned DIM>
362 {
363  return mNumBoxes;
364 }
365 
366 template<unsigned DIM>
368 {
369  return mBoxes.size();
370 }
371 
372 template<unsigned DIM>
373 c_vector<double, 2*DIM> DistributedBoxCollection<DIM>::rGetDomainSize() const
374 {
375  return mDomainSize;
376 }
377 
378 template<unsigned DIM>
380 {
381  return mAreLocalBoxesSet;
382 }
383 
384 template<unsigned DIM>
386 {
387  return mBoxWidth;
388 }
389 
390 template<unsigned DIM>
392 {
393  return mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow();
394 }
395 
396 template<unsigned DIM>
397 int DistributedBoxCollection<DIM>::LoadBalance(std::vector<int> localDistribution)
398 {
399  MPI_Status status;
400 
401  int proc_right = (PetscTools::AmTopMost()) ? MPI_PROC_NULL : (int)PetscTools::GetMyRank() + 1;
402  int proc_left = (PetscTools::AmMaster()) ? MPI_PROC_NULL : (int)PetscTools::GetMyRank() - 1;
403 
404  // A variable that will return the new number of rows.
405  int new_rows = localDistribution.size();
406 
410  unsigned rows_on_left_process = 0;
411  std::vector<int> node_distr_on_left_process;
412 
413  unsigned num_local_rows = localDistribution.size();
414 
415  MPI_Send(&num_local_rows, 1, MPI_UNSIGNED, proc_right, 123, PETSC_COMM_WORLD);
416  MPI_Recv(&rows_on_left_process, 1, MPI_UNSIGNED, proc_left, 123, PETSC_COMM_WORLD, &status);
417 
418  node_distr_on_left_process.resize(rows_on_left_process > 0 ? rows_on_left_process : 1);
419 
420  MPI_Send(&localDistribution[0], num_local_rows, MPI_INT, proc_right, 123, PETSC_COMM_WORLD);
421  MPI_Recv(&node_distr_on_left_process[0], rows_on_left_process, MPI_INT, proc_left, 123, PETSC_COMM_WORLD, &status);
422 
426  int local_load = 0;
427  for (unsigned i=0; i<localDistribution.size(); i++)
428  {
429  local_load += localDistribution[i];
430  }
431  int load_on_left_proc = 0;
432  for (unsigned i=0; i<node_distr_on_left_process.size(); i++)
433  {
434  load_on_left_proc += node_distr_on_left_process[i];
435  }
436 
437  if (!PetscTools::AmMaster())
438  {
439  // Calculate (Difference in load with a shift) - (Difference in current loads) for a move left and right of the boundary
440  // This code uses integer arithmetic in order to avoid the rounding errors associated with doubles
441  int local_to_left_sq = (local_load - load_on_left_proc) * (local_load - load_on_left_proc);
442  int delta_left = ( (local_load + node_distr_on_left_process[node_distr_on_left_process.size() - 1]) - (load_on_left_proc - node_distr_on_left_process[node_distr_on_left_process.size() - 1]) );
443  delta_left = delta_left*delta_left - local_to_left_sq;
444 
445  int delta_right = ( (local_load - localDistribution[0]) - (load_on_left_proc + localDistribution[0]));
446  delta_right = delta_right*delta_right - local_to_left_sq;
447 
448  // If a delta is negative we should accept that change. If both are negative choose the largest change.
449  int local_change = 0;
450  bool move_left = (!(delta_left > 0) && (node_distr_on_left_process.size() > 1));
451  if (move_left)
452  {
453  local_change = 1;
454  }
455 
456  bool move_right = !(delta_right > 0) && (localDistribution.size() > 2);
457  if (move_right)
458  {
459  local_change = -1;
460  }
461 
462  if (move_left && move_right)
463  {
464  local_change = (fabs((double)delta_right) > fabs((double)delta_left)) ? -1 : 1;
465  }
466 
467  // Update the number of local rows.
468  new_rows += local_change;
469 
470  // Send the result of the calculation back to the left processes.
471  MPI_Send(&local_change, 1, MPI_INT, proc_left, 123, PETSC_COMM_WORLD);
472  }
473 
474  // Receive changes from right hand process.
475  int remote_change = 0;
476  MPI_Recv(&remote_change, 1, MPI_INT, proc_right, 123, PETSC_COMM_WORLD, &status);
477 
478  // Update based on change or right/top boundary
479  new_rows -= remote_change;
480 
481  return new_rows;
482 }
483 
484 template<unsigned DIM>
486 {
487  if (mAreLocalBoxesSet)
488  {
489  EXCEPTION("Local Boxes Are Already Set");
490  }
491  else
492  {
493  switch (DIM)
494  {
495  case 1:
496  {
497  // We only need to look for neighbours in the current and successive boxes plus some others for halos
498  mLocalBoxes.clear();
499 
500  // Iterate over the global box indices
501  for (unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
502  {
503  std::set<unsigned> local_boxes;
504 
505  // Insert the current box
506  local_boxes.insert(global_index);
507 
508  // Set some bools to find out where we are
509  bool right = (global_index==mNumBoxesEachDirection(0)-1);
510  bool left = (global_index == 0);
511  bool proc_left = (global_index == mpDistributedBoxStackFactory->GetLow());
512 
513  // If we're not at the right-most box, then insert the box to the right
514  if (!right)
515  {
516  local_boxes.insert(global_index+1);
517  }
518  // If we're on a left process boundary and not on process 0, insert the (halo) box to the left
519  if (proc_left && !left)
520  {
521  local_boxes.insert(global_index-1);
522  }
523 
524  mLocalBoxes.push_back(local_boxes);
525  }
526  break;
527  }
528  case 2:
529  {
530  // We only need to look for neighbours in the current box and half the neighbouring boxes plus some others for halos
531  mLocalBoxes.clear();
532 
533  for (unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
534  {
535  std::set<unsigned> local_boxes;
536 
537  // Set up bools to find out where we are
538  bool left = (global_index%mNumBoxesEachDirection(0) == 0);
539  bool right = (global_index%mNumBoxesEachDirection(0) == mNumBoxesEachDirection(0)-1);
540  bool top = !(global_index < mNumBoxesEachDirection(0)*mNumBoxesEachDirection(1) - mNumBoxesEachDirection(0));
541  bool bottom = (global_index < mNumBoxesEachDirection(0));
542  bool bottom_proc = (CalculateCoordinateIndices(global_index)[1] == mpDistributedBoxStackFactory->GetLow());
543 
544  // Insert the current box
545  local_boxes.insert(global_index);
546 
547  // If we're on the bottom of the process boundary, but not the bottom of the domain add boxes below
548  if(!bottom && bottom_proc)
549  {
550  local_boxes.insert(global_index - mNumBoxesEachDirection(0));
551  if(!left)
552  {
553  local_boxes.insert(global_index - mNumBoxesEachDirection(0) - 1);
554  }
555  if(!right)
556  {
557  local_boxes.insert(global_index - mNumBoxesEachDirection(0) + 1);
558  }
559  }
560 
561  // If we're not at the top of the domain insert boxes above
562  if(!top)
563  {
564  local_boxes.insert(global_index + mNumBoxesEachDirection(0));
565 
566  if(!right)
567  {
568  local_boxes.insert(global_index + mNumBoxesEachDirection(0) + 1);
569  }
570  if(!left)
571  {
572  local_boxes.insert(global_index + mNumBoxesEachDirection(0) - 1);
573  }
574  else if ( (global_index % mNumBoxesEachDirection(0) == 0) && (mIsPeriodicInX) ) // If we're on the left edge but its periodic include the box on the far right and up one.
575  {
576  local_boxes.insert(global_index + 2 * mNumBoxesEachDirection(0) - 1);
577  }
578  }
579 
580  // If we're not on the far right hand side inseryt box to the right
581  if(!right)
582  {
583  local_boxes.insert(global_index + 1);
584  }
585  // If we're on the right edge but it's periodic include the box on the far left of the domain
586  else if ( (global_index % mNumBoxesEachDirection(0) == mNumBoxesEachDirection(0)-1) && (mIsPeriodicInX) )
587  {
588  local_boxes.insert(global_index - mNumBoxesEachDirection(0) + 1);
589  // If we're also not on the top-most row, then insert the box above- on the far left of the domain
590  if (global_index < mBoxes.size() - mNumBoxesEachDirection(0))
591  {
592  local_boxes.insert(global_index + 1);
593  }
594  }
595 
596  mLocalBoxes.push_back(local_boxes);
597  }
598  break;
599  }
600  case 3:
601  {
602  // We only need to look for neighbours in the current box and half the neighbouring boxes plus some others for halos
603  mLocalBoxes.clear();
604  unsigned num_boxes_xy = mNumBoxesEachDirection(0)*mNumBoxesEachDirection(1);
605 
606 
607  for (unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
608  {
609  std::set<unsigned> local_boxes;
610 
611  // Set some bools to find out where we are
612  bool top = !(global_index % num_boxes_xy < num_boxes_xy - mNumBoxesEachDirection(0));
613  bool bottom = (global_index % num_boxes_xy < mNumBoxesEachDirection(0));
614  bool left = (global_index % mNumBoxesEachDirection(0) == 0);
615  bool right = (global_index % mNumBoxesEachDirection(0) == mNumBoxesEachDirection(0) - 1);
616  bool front = (global_index < num_boxes_xy);
617  bool back = !(global_index < num_boxes_xy*mNumBoxesEachDirection(2) - num_boxes_xy);
618  bool proc_front = (CalculateCoordinateIndices(global_index)[2] == mpDistributedBoxStackFactory->GetLow());
619  bool proc_back = (CalculateCoordinateIndices(global_index)[2] == mpDistributedBoxStackFactory->GetHigh()-1);
620 
621  // Insert the current box
622  local_boxes.insert(global_index);
623 
624  // If we're not on the front face, add appropriate boxes on the closer face
625  if(!front)
626  {
627  // If we're not on the top of the domain
628  if(!top)
629  {
630  local_boxes.insert( global_index - num_boxes_xy + mNumBoxesEachDirection(0) );
631  if(!left)
632  {
633  local_boxes.insert( global_index - num_boxes_xy + mNumBoxesEachDirection(0) - 1);
634  }
635  if(!right)
636  {
637  local_boxes.insert( global_index - num_boxes_xy + mNumBoxesEachDirection(0) + 1);
638  }
639  }
640  if(!right)
641  {
642  local_boxes.insert( global_index - num_boxes_xy + 1);
643  }
644 
645  // If we are on the front of the process we have to add extra boxes as they are halos.
646  if(proc_front)
647  {
648  local_boxes.insert( global_index - num_boxes_xy );
649 
650  if(!left)
651  {
652  local_boxes.insert( global_index - num_boxes_xy - 1);
653  }
654  if(!bottom)
655  {
656  local_boxes.insert( global_index - num_boxes_xy - mNumBoxesEachDirection(0));
657 
658  if(!left)
659  {
660  local_boxes.insert( global_index - num_boxes_xy - mNumBoxesEachDirection(0) - 1);
661  }
662  if(!right)
663  {
664  local_boxes.insert( global_index - num_boxes_xy - mNumBoxesEachDirection(0) + 1);
665  }
666  }
667 
668  }
669  }
670  if(!right)
671  {
672  local_boxes.insert( global_index + 1);
673  }
674  // If we're not on the very top add boxes above
675  if(!top)
676  {
677  local_boxes.insert( global_index + mNumBoxesEachDirection(0));
678 
679  if(!right)
680  {
681  local_boxes.insert( global_index + mNumBoxesEachDirection(0) + 1);
682  }
683  if(!left)
684  {
685  local_boxes.insert( global_index + mNumBoxesEachDirection(0) - 1);
686  }
687  }
688 
689  // If we're not on the back add boxes behind
690  if(!back)
691  {
692  local_boxes.insert(global_index + num_boxes_xy);
693 
694  if(!right)
695  {
696  local_boxes.insert(global_index + num_boxes_xy + 1);
697  }
698  if(!top)
699  {
700  local_boxes.insert(global_index + num_boxes_xy + mNumBoxesEachDirection(0));
701  if(!right)
702  {
703  local_boxes.insert(global_index + num_boxes_xy + mNumBoxesEachDirection(0) + 1);
704  }
705  if(!left)
706  {
707  local_boxes.insert(global_index + num_boxes_xy + mNumBoxesEachDirection(0) - 1);
708  }
709  }
710  // If we are on the back proc we should make sure we get everything in the face further back
711  if(proc_back)
712  {
713  if(!left)
714  {
715  local_boxes.insert(global_index + num_boxes_xy - 1);
716  }
717  if(!bottom)
718  {
719  local_boxes.insert(global_index + num_boxes_xy - mNumBoxesEachDirection(0));
720 
721  if(!left)
722  {
723  local_boxes.insert(global_index + num_boxes_xy - mNumBoxesEachDirection(0) - 1);
724  }
725  if(!right)
726  {
727  local_boxes.insert(global_index + num_boxes_xy - mNumBoxesEachDirection(0) + 1);
728  }
729  }
730  }
731  }
732  mLocalBoxes.push_back(local_boxes);
733  }
734 
735  break;
736  }
737  default:
739  }
740  mAreLocalBoxesSet=true;
741  }
742 }
743 
744 
745 
746 template<unsigned DIM>
748 {
749  mAreLocalBoxesSet = true;
750  switch (DIM)
751  {
752  case 1:
753  {
754  for (unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
755  {
756  std::set<unsigned> local_boxes;
757 
758  local_boxes.insert(i);
759 
760  // add the two neighbours
761  if (i!=0)
762  {
763  local_boxes.insert(i-1);
764  }
765  if (i+1 != mNumBoxesEachDirection(0))
766  {
767  local_boxes.insert(i+1);
768  }
769 
770  mLocalBoxes.push_back(local_boxes);
771  }
772  break;
773  }
774  case 2:
775  {
776  mLocalBoxes.clear();
777 
778  unsigned M = mNumBoxesEachDirection(0);
779  unsigned N = mNumBoxesEachDirection(1);
780 
781  std::vector<bool> is_xmin(N*M); // far left
782  std::vector<bool> is_xmax(N*M); // far right
783  std::vector<bool> is_ymin(N*M); // bottom
784  std::vector<bool> is_ymax(N*M); // top
785 
786  for (unsigned i=0; i<M*N; i++)
787  {
788  is_xmin[i] = (i%M==0);
789  is_xmax[i] = ((i+1)%M==0);
790  is_ymin[i] = (i%(M*N)<M);
791  is_ymax[i] = (i%(M*N)>=(N-1)*M);
792  }
793 
794  for (unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
795  {
796  std::set<unsigned> local_boxes;
797 
798  local_boxes.insert(i);
799 
800  // add the box to the left
801  if (!is_xmin[i])
802  {
803  local_boxes.insert(i-1);
804  }
805  else // Add Periodic Box if needed
806  {
807  if(mIsPeriodicInX)
808  {
809  local_boxes.insert(i+M-1);
810  }
811  }
812 
813  // add the box to the right
814  if (!is_xmax[i])
815  {
816  local_boxes.insert(i+1);
817  }
818  else // Add Periodic Box if needed
819  {
820  if(mIsPeriodicInX)
821  {
822  local_boxes.insert(i-M+1);
823  }
824  }
825 
826  // add the one below
827  if (!is_ymin[i])
828  {
829  local_boxes.insert(i-M);
830  }
831 
832  // add the one above
833  if (!is_ymax[i])
834  {
835  local_boxes.insert(i+M);
836  }
837 
838  // add the four corner boxes
839 
840  if ( (!is_xmin[i]) && (!is_ymin[i]) )
841  {
842  local_boxes.insert(i-1-M);
843  }
844  if ( (!is_xmin[i]) && (!is_ymax[i]) )
845  {
846  local_boxes.insert(i-1+M);
847  }
848  if ( (!is_xmax[i]) && (!is_ymin[i]) )
849  {
850  local_boxes.insert(i+1-M);
851  }
852  if ( (!is_xmax[i]) && (!is_ymax[i]) )
853  {
854  local_boxes.insert(i+1+M);
855  }
856 
857  // Add Periodic Corner Boxes if needed
858  if(mIsPeriodicInX)
859  {
860  if( (is_xmin[i]) && (!is_ymin[i]) )
861  {
862  local_boxes.insert(i-1);
863  }
864  if ( (is_xmin[i]) && (!is_ymax[i]) )
865  {
866  local_boxes.insert(i-1+2*M);
867  }
868  if ( (is_xmax[i]) && (!is_ymin[i]) )
869  {
870  local_boxes.insert(i+1-2*M);
871  }
872  if ( (is_xmax[i]) && (!is_ymax[i]) )
873  {
874  local_boxes.insert(i+1);
875  }
876  }
877 
878  mLocalBoxes.push_back(local_boxes);
879  }
880  break;
881  }
882  case 3:
883  {
884  mLocalBoxes.clear();
885 
886  unsigned M = mNumBoxesEachDirection(0);
887  unsigned N = mNumBoxesEachDirection(1);
888  unsigned P = mNumBoxesEachDirection(2);
889 
890  std::vector<bool> is_xmin(N*M*P); // far left
891  std::vector<bool> is_xmax(N*M*P); // far right
892  std::vector<bool> is_ymin(N*M*P); // nearest
893  std::vector<bool> is_ymax(N*M*P); // furthest
894  std::vector<bool> is_zmin(N*M*P); // bottom layer
895  std::vector<bool> is_zmax(N*M*P); // top layer
896 
897  for (unsigned i=0; i<M*N*P; i++)
898  {
899  is_xmin[i] = (i%M==0);
900  is_xmax[i] = ((i+1)%M==0);
901  is_ymin[i] = (i%(M*N)<M);
902  is_ymax[i] = (i%(M*N)>=(N-1)*M);
903  is_zmin[i] = (i<M*N);
904  is_zmax[i] = (i>=M*N*(P-1));
905  }
906 
907  for (unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
908  {
909  std::set<unsigned> local_boxes;
910 
911  // add itself as a local box
912  local_boxes.insert(i);
913 
914  // now add all 26 other neighbours.....
915 
916  // add the box left
917  if (!is_xmin[i])
918  {
919  local_boxes.insert(i-1);
920 
921  // plus some others towards the left
922  if (!is_ymin[i])
923  {
924  local_boxes.insert(i-1-M);
925  }
926 
927  if (!is_ymax[i])
928  {
929  local_boxes.insert(i-1+M);
930  }
931 
932  if (!is_zmin[i])
933  {
934  local_boxes.insert(i-1-M*N);
935  }
936 
937  if (!is_zmax[i])
938  {
939  local_boxes.insert(i-1+M*N);
940  }
941  }
942 
943  // add the box to the right
944  if (!is_xmax[i])
945  {
946  local_boxes.insert(i+1);
947 
948  // plus some others towards the right
949  if (!is_ymin[i])
950  {
951  local_boxes.insert(i+1-M);
952  }
953 
954  if (!is_ymax[i])
955  {
956  local_boxes.insert(i+1+M);
957  }
958 
959  if (!is_zmin[i])
960  {
961  local_boxes.insert(i+1-M*N);
962  }
963 
964  if (!is_zmax[i])
965  {
966  local_boxes.insert(i+1+M*N);
967  }
968  }
969 
970  // add the boxes next along the y axis
971  if (!is_ymin[i])
972  {
973  local_boxes.insert(i-M);
974 
975  // and more in this plane
976  if (!is_zmin[i])
977  {
978  local_boxes.insert(i-M-M*N);
979  }
980 
981  if (!is_zmax[i])
982  {
983  local_boxes.insert(i-M+M*N);
984  }
985  }
986 
987  // add the boxes next along the y axis
988  if (!is_ymax[i])
989  {
990  local_boxes.insert(i+M);
991 
992  // and more in this plane
993  if (!is_zmin[i])
994  {
995  local_boxes.insert(i+M-M*N);
996  }
997 
998  if (!is_zmax[i])
999  {
1000  local_boxes.insert(i+M+M*N);
1001  }
1002  }
1003 
1004  // add the box directly above
1005  if (!is_zmin[i])
1006  {
1007  local_boxes.insert(i-N*M);
1008  }
1009 
1010  // add the box directly below
1011  if (!is_zmax[i])
1012  {
1013  local_boxes.insert(i+N*M);
1014  }
1015 
1016  // finally, the 8 corners are left
1017 
1018  if ( (!is_xmin[i]) && (!is_ymin[i]) && (!is_zmin[i]) )
1019  {
1020  local_boxes.insert(i-1-M-M*N);
1021  }
1022 
1023  if ( (!is_xmin[i]) && (!is_ymin[i]) && (!is_zmax[i]) )
1024  {
1025  local_boxes.insert(i-1-M+M*N);
1026  }
1027 
1028  if ( (!is_xmin[i]) && (!is_ymax[i]) && (!is_zmin[i]) )
1029  {
1030  local_boxes.insert(i-1+M-M*N);
1031  }
1032 
1033  if ( (!is_xmin[i]) && (!is_ymax[i]) && (!is_zmax[i]) )
1034  {
1035  local_boxes.insert(i-1+M+M*N);
1036  }
1037 
1038  if ( (!is_xmax[i]) && (!is_ymin[i]) && (!is_zmin[i]) )
1039  {
1040  local_boxes.insert(i+1-M-M*N);
1041  }
1042 
1043  if ( (!is_xmax[i]) && (!is_ymin[i]) && (!is_zmax[i]) )
1044  {
1045  local_boxes.insert(i+1-M+M*N);
1046  }
1047 
1048  if ( (!is_xmax[i]) && (!is_ymax[i]) && (!is_zmin[i]) )
1049  {
1050  local_boxes.insert(i+1+M-M*N);
1051  }
1052 
1053  if ( (!is_xmax[i]) && (!is_ymax[i]) && (!is_zmax[i]) )
1054  {
1055  local_boxes.insert(i+1+M+M*N);
1056  }
1057 
1058  mLocalBoxes.push_back(local_boxes);
1059  }
1060  break;
1061  }
1062  default:
1063  NEVER_REACHED;
1064  }
1065 }
1066 
1067 template<unsigned DIM>
1068 std::set<unsigned> DistributedBoxCollection<DIM>::GetLocalBoxes(unsigned boxIndex)
1069 {
1070  // Make sure the box is locally owned
1071  assert(!(boxIndex < mMinBoxIndex) && !(mMaxBoxIndex<boxIndex));
1072  return mLocalBoxes[boxIndex-mMinBoxIndex];
1073 }
1074 
1075 template<unsigned DIM>
1077 {
1078  unsigned index = CalculateContainingBox(pNode);
1079 
1080  return GetBoxOwnership(index);
1081 }
1082 
1083 template<unsigned DIM>
1084 bool DistributedBoxCollection<DIM>::IsOwned(c_vector<double, DIM>& location)
1085 {
1086  unsigned index = CalculateContainingBox(location);
1087 
1088  return GetBoxOwnership(index);
1089 }
1090 
1091 template<unsigned DIM>
1093 {
1094  unsigned box_index = CalculateContainingBox(pNode);
1095  unsigned containing_process = PetscTools::GetMyRank();
1096 
1097  if (box_index > mMaxBoxIndex)
1098  {
1099  containing_process++;
1100  }
1101  else if (box_index < mMinBoxIndex)
1102  {
1103  containing_process--;
1104  }
1105 
1106  return containing_process;
1107 }
1108 
1109 template<unsigned DIM>
1111 {
1112  return mHaloNodesRight;
1113 }
1114 
1115 template<unsigned DIM>
1117 {
1118  return mHaloNodesLeft;
1119 }
1120 
1121 template<unsigned DIM>
1123 {
1124  mCalculateNodeNeighbours = calculateNodeNeighbours;
1125 }
1126 
1127 template<unsigned DIM>
1128 void DistributedBoxCollection<DIM>::CalculateNodePairs(std::vector<Node<DIM>*>& rNodes, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs, std::map<unsigned, std::set<unsigned> >& rNodeNeighbours)
1129 {
1130  rNodePairs.clear();
1131  rNodeNeighbours.clear();
1132 
1133  // Create an empty neighbours set for each node
1134  for (unsigned i=0; i<rNodes.size(); i++)
1135  {
1136  unsigned node_index = rNodes[i]->GetIndex();
1137 
1138  // Get the box containing this node
1139  unsigned box_index = CalculateContainingBox(rNodes[i]);
1140 
1141  if (GetBoxOwnership(box_index))
1142  {
1143  rNodeNeighbours[node_index] = std::set<unsigned>();
1144  }
1145  }
1146 
1147  for (std::map<unsigned, unsigned>::iterator map_iter = mBoxesMapping.begin();
1148  map_iter != mBoxesMapping.end();
1149  ++map_iter)
1150  {
1151  // Get the box global index
1152  unsigned box_index = map_iter->first;
1153 
1154  AddPairsFromBox(box_index, rNodePairs, rNodeNeighbours);
1155  }
1156 }
1157 
1158 template<unsigned DIM>
1159 void DistributedBoxCollection<DIM>::CalculateInteriorNodePairs(std::vector<Node<DIM>*>& rNodes, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs, std::map<unsigned, std::set<unsigned> >& rNodeNeighbours)
1160 {
1161  rNodePairs.clear();
1162  rNodeNeighbours.clear();
1163 
1164  // Create an empty neighbours set for each node
1165  for (unsigned i=0; i<rNodes.size(); i++)
1166  {
1167  unsigned node_index = rNodes[i]->GetIndex();
1168 
1169  // Get the box containing this node
1170  unsigned box_index = CalculateContainingBox(rNodes[i]);
1171 
1172  if (GetBoxOwnership(box_index))
1173  {
1174  rNodeNeighbours[node_index] = std::set<unsigned>();
1175  }
1176  }
1177 
1178  for (std::map<unsigned, unsigned>::iterator map_iter = mBoxesMapping.begin();
1179  map_iter != mBoxesMapping.end();
1180  ++map_iter)
1181  {
1182  // Get the box global index
1183  unsigned box_index = map_iter->first;
1184 
1185  if (IsInteriorBox(box_index))
1186  {
1187  AddPairsFromBox(box_index, rNodePairs, rNodeNeighbours);
1188  }
1189  }
1190 }
1191 
1192 template<unsigned DIM>
1193 void DistributedBoxCollection<DIM>::CalculateBoundaryNodePairs(std::vector<Node<DIM>*>& rNodes, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs, std::map<unsigned, std::set<unsigned> >& rNodeNeighbours)
1194 {
1195  for (std::map<unsigned, unsigned>::iterator map_iter = mBoxesMapping.begin();
1196  map_iter != mBoxesMapping.end();
1197  ++map_iter)
1198  {
1199  // Get the box global index
1200  unsigned box_index = map_iter->first;
1201 
1202  if (!IsInteriorBox(box_index))
1203  {
1204  AddPairsFromBox(box_index, rNodePairs, rNodeNeighbours);
1205  }
1206  }
1207 }
1208 
1209 template<unsigned DIM>
1210 void DistributedBoxCollection<DIM>::AddPairsFromBox(unsigned boxIndex, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs, std::map<unsigned, std::set<unsigned> >& rNodeNeighbours)
1211 {
1212  // Get the box
1213  Box<DIM> box = rGetBox(boxIndex);
1214 
1215  // Get the set of nodes in this box
1216  std::set< Node<DIM>* >& r_contained_nodes = box.rGetNodesContained();
1217 
1218  // Get the local boxes to this box
1219  std::set<unsigned> local_boxes_indices = GetLocalBoxes(boxIndex);
1220 
1221  // Loop over all the local boxes
1222  for (std::set<unsigned>::iterator box_iter = local_boxes_indices.begin();
1223  box_iter != local_boxes_indices.end();
1224  box_iter++)
1225  {
1226  Box<DIM>* p_neighbour_box;
1227 
1228  // Establish whether box is locally owned or halo.
1229  if (GetBoxOwnership(*box_iter))
1230  {
1231  p_neighbour_box = &mBoxes[mBoxesMapping[*box_iter]];
1232  }
1233  else // Assume it is a halo.
1234  {
1235  p_neighbour_box = &mHaloBoxes[mHaloBoxesMapping[*box_iter]];
1236  }
1237  assert(p_neighbour_box);
1238 
1239  // Get the set of nodes contained in this box
1240  std::set< Node<DIM>* >& r_contained_neighbour_nodes = p_neighbour_box->rGetNodesContained();
1241 
1242  // Loop over these nodes
1243  for (typename std::set<Node<DIM>*>::iterator neighbour_node_iter = r_contained_neighbour_nodes.begin();
1244  neighbour_node_iter != r_contained_neighbour_nodes.end();
1245  ++neighbour_node_iter)
1246  {
1247  // Get the index of the other node
1248  unsigned other_node_index = (*neighbour_node_iter)->GetIndex();
1249 
1250  // Loop over nodes in this box
1251  for (typename std::set<Node<DIM>*>::iterator node_iter = r_contained_nodes.begin();
1252  node_iter != r_contained_nodes.end();
1253  ++node_iter)
1254  {
1255  unsigned node_index = (*node_iter)->GetIndex();
1256 
1257  // If we're in the same box, then take care not to store the node pair twice
1258  if (*box_iter == boxIndex)
1259  {
1260  if (other_node_index > node_index)
1261  {
1262  rNodePairs.push_back(std::pair<Node<DIM>*, Node<DIM>*>((*node_iter), (*neighbour_node_iter)));
1263  if (mCalculateNodeNeighbours)
1264  {
1265  rNodeNeighbours[node_index].insert(other_node_index);
1266  rNodeNeighbours[other_node_index].insert(node_index);
1267  }
1268  }
1269  }
1270  else
1271  {
1272  rNodePairs.push_back(std::pair<Node<DIM>*, Node<DIM>*>((*node_iter), (*neighbour_node_iter)));
1273  if (mCalculateNodeNeighbours)
1274  {
1275  rNodeNeighbours[node_index].insert(other_node_index);
1276  rNodeNeighbours[other_node_index].insert(node_index);
1277  }
1278  }
1279 
1280  }
1281  }
1282  }
1283 }
1284 
1285 template<unsigned DIM>
1287 {
1288  std::vector<int> cell_numbers(mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow(), 0);
1289 
1290  for (std::map<unsigned, unsigned>::iterator iter = mBoxesMapping.begin();
1291  iter != mBoxesMapping.end();
1292  ++iter)
1293  {
1294  c_vector<unsigned, DIM> coords = CalculateCoordinateIndices(iter->first);
1295  unsigned location_in_vector = coords[DIM-1]-mpDistributedBoxStackFactory->GetLow();
1296 
1297  cell_numbers[location_in_vector] += mBoxes[iter->second].rGetNodesContained().size();
1298  }
1299 
1300  return cell_numbers;
1301 }
1302 
1304 // Explicit instantiation
1306 
1307 template class DistributedBoxCollection<1>;
1308 template class DistributedBoxCollection<2>;
1309 template class DistributedBoxCollection<3>;
std::vector< Box< DIM > > mBoxes
std::vector< unsigned > & rGetHaloNodesRight()
std::set< unsigned > GetLocalBoxes(unsigned boxIndex)
unsigned CalculateContainingBox(Node< DIM > *pNode)
void CalculateInteriorNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs, std::map< unsigned, std::set< unsigned > > &rNodeNeighbours)
Definition: Node.hpp:58
bool IsInteriorBox(unsigned globalIndex)
c_vector< double, 2 *DIM > rGetDomainSize() const
c_vector< double, 2 *DIM > mDomainSize
void AddPairsFromBox(unsigned boxIndex, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs, std::map< unsigned, std::set< unsigned > > &rNodeNeighbours)
#define EXCEPTION(message)
Definition: Exception.hpp:143
void SetCalculateNodeNeighbours(bool calculateNodeNeighbours)
static bool AmMaster()
Definition: PetscTools.cpp:120
std::set< Node< DIM > * > & rGetNodesContained()
Definition: Box.cpp:73
std::map< unsigned, unsigned > mBoxesMapping
std::vector< int > CalculateNumberOfNodesInEachStrip()
unsigned CalculateGlobalIndex(c_vector< unsigned, DIM > coordinateIndices)
#define NEVER_REACHED
Definition: Exception.hpp:206
unsigned GetProcessOwningNode(Node< DIM > *pNode)
static bool IsSequential()
Definition: PetscTools.cpp:91
c_vector< unsigned, DIM > CalculateCoordinateIndices(unsigned globalIndex)
c_vector< unsigned, DIM > mNumBoxesEachDirection
std::vector< unsigned > & rGetHaloNodesLeft()
bool IsOwned(Node< DIM > *pNode)
DistributedBoxCollection(double boxWidth, c_vector< double, 2 *DIM > domainSize, bool isPeriodicInX=false, int localRows=PETSC_DECIDE)
bool GetHaloBoxOwnership(unsigned globalIndex)
bool GetBoxOwnership(unsigned globalIndex)
void CalculateNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs, std::map< unsigned, std::set< unsigned > > &rNodeNeighbours)
int LoadBalance(std::vector< int > localDistribution)
Definition: Box.hpp:50
static bool IsParallel()
Definition: PetscTools.cpp:97
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition: Node.cpp:140
Box< DIM > & rGetHaloBox(unsigned boxIndex)
static bool AmTopMost()
Definition: PetscTools.cpp:126
static unsigned GetNumProcs()
Definition: PetscTools.cpp:108
static unsigned GetMyRank()
Definition: PetscTools.cpp:114
Box< DIM > & rGetBox(unsigned boxIndex)
DistributedVectorFactory * mpDistributedBoxStackFactory
void CalculateBoundaryNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs, std::map< unsigned, std::set< unsigned > > &rNodeNeighbours)