Chaste Commit::f2ff7ee04e70ac9d06c57344df8d017dbb12b97b
DistributedBoxCollection.cpp
1/*
2
3Copyright (c) 2005-2024, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, 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
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35#include "DistributedBoxCollection.hpp"
36#include "Exception.hpp"
38#include "Warnings.hpp"
39
40// Static member for "fudge factor" is instantiated here
41template<unsigned DIM>
43
44template<unsigned DIM>
45DistributedBoxCollection<DIM>::DistributedBoxCollection(double boxWidth, c_vector<double, 2*DIM> domainSize, bool isPeriodicInX, bool isPeriodicInY, bool isPeriodicInZ, int localRows)
46 : mBoxWidth(boxWidth),
47 mIsPeriodicInX(isPeriodicInX),
48 mIsPeriodicInY(isPeriodicInY),
49 mIsPeriodicInZ(isPeriodicInZ),
50 mAreLocalBoxesSet(false),
51 mCalculateNodeNeighbours(true)
52{
53 // If the domain size is not 'divisible' (i.e. fmod(width, box_size) > 0.0) we swell the domain to enforce this.
54 for (unsigned i=0; i<DIM; i++)
55 {
56 double r = fmod((domainSize[2*i+1]-domainSize[2*i]), boxWidth);
57 if (r > 0.0)
58 {
59 domainSize[2*i+1] += boxWidth - r;
60 }
61 }
62
63 mDomainSize = domainSize;
64
65 // Set the periodicity across procs flag
66 mIsPeriodicAcrossProcs = (DIM==1 && mIsPeriodicInX) || (DIM==2 && mIsPeriodicInY) || (DIM==3 && mIsPeriodicInZ);
67
68 // Calculate the number of boxes in each direction.
69 mNumBoxesEachDirection = scalar_vector<unsigned>(DIM, 0u);
70
71 for (unsigned i=0; i<DIM; i++)
72 {
73 double counter = mDomainSize(2*i);
74 while (counter + msFudge < mDomainSize(2*i+1))
75 {
77 counter += mBoxWidth;
78 }
79 }
80
81 // Make sure there are enough boxes for the number of processes.
83 {
84 WARNING("There are more processes than convenient for the domain/mesh/box size. The domain size has been swollen.")
87 }
88
89 // Make a distributed vector factory to split the rows of boxes between processes.
91
92 // Calculate how many boxes in a row / face. A useful piece of data in the class.
93 mNumBoxes = 1u;
94 for (unsigned dim=0; dim<DIM; dim++)
95 {
97 }
98
100
101 unsigned num_local_boxes = mNumBoxesInAFace * GetNumLocalRows();
102
105
106 // Create the correct number of boxes and set up halos
107 mBoxes.resize(num_local_boxes);
109}
110
111template<unsigned DIM>
113{
114 delete mpDistributedBoxStackFactory;
115}
116
117template<unsigned DIM>
119{
120 for (unsigned i=0; i<mBoxes.size(); i++)
121 {
122 mBoxes[i].ClearNodes();
123 }
124 for (unsigned i=0; i<mHaloBoxes.size(); i++)
125 {
126 mHaloBoxes[i].ClearNodes();
127 }
128}
129
130template<unsigned DIM>
132{
133 // We don't need to do this if not parallel
135 {
136 return;
137 }
138
139 // Get top-most and bottom-most value of Distributed Box Stack.
140 unsigned hi = mpDistributedBoxStackFactory->GetHigh();
141 unsigned lo = mpDistributedBoxStackFactory->GetLow();
142
143 // If I am not the top-most process, add halo structures above.
145 {
146 for (unsigned i=0; i < mNumBoxesInAFace; i++)
147 {
148 Box<DIM> new_box;
149 mHaloBoxes.push_back(new_box);
150
151 unsigned global_index = hi * mNumBoxesInAFace + i;
152 mHaloBoxesMapping[global_index] = mHaloBoxes.size()-1;
153 mHalosRight.push_back(global_index - mNumBoxesInAFace);
154 }
155 }
156 // Otherwise if I am the top most and periodic in y (2d) or z (3d) add halo boxes for
157 // the base process
158 else if ( mIsPeriodicAcrossProcs )
159 {
160 for (unsigned i=0; i < mNumBoxesInAFace; i++)
161 {
162 Box<DIM> new_box;
163 mHaloBoxes.push_back(new_box);
164
165 mHaloBoxesMapping[i] = mHaloBoxes.size()-1;
167 mHalosRight.push_back( (hi-1)*mNumBoxesInAFace + i );
168 }
169 }
170
171 // If I am not the bottom-most process, add halo structures below.
173 {
174 for (unsigned i=0; i< mNumBoxesInAFace; i++)
175 {
176 Box<DIM> new_box;
177 mHaloBoxes.push_back(new_box);
178
179 unsigned global_index = (lo - 1) * mNumBoxesInAFace + i;
180 mHaloBoxesMapping[global_index] = mHaloBoxes.size() - 1;
181
182 mHalosLeft.push_back(global_index + mNumBoxesInAFace);
183 }
184 }
185 // Otherwise if I am the bottom most and periodic in y (2d) or z (3d) add halo boxes for
186 // the top process
187 else if ( mIsPeriodicAcrossProcs )
188 {
189 for (unsigned i=0; i < mNumBoxesInAFace; i++)
190 {
191 Box<DIM> new_box;
192 mHaloBoxes.push_back(new_box);
194 unsigned global_index = (mNumBoxesEachDirection(DIM-1) - 1) * mNumBoxesInAFace + i;
195 mHaloBoxesMapping[global_index] = mHaloBoxes.size()-1;
196
197 mHalosLeft.push_back( i );
198 }
199
200 }
201}
203template<unsigned DIM>
205{
206 mHaloNodesLeft.clear();
207 for (unsigned i=0; i<mHalosLeft.size(); i++)
209 for (typename std::set<Node<DIM>* >::iterator iter=this->rGetBox(mHalosLeft[i]).rGetNodesContained().begin();
210 iter!=this->rGetBox(mHalosLeft[i]).rGetNodesContained().end();
211 iter++)
212 {
213 mHaloNodesLeft.push_back((*iter)->GetIndex());
214 }
215 }
217 // Send right
218 mHaloNodesRight.clear();
219 for (unsigned i=0; i<mHalosRight.size(); i++)
220 {
221 for (typename std::set<Node<DIM>* >::iterator iter=this->rGetBox(mHalosRight[i]).rGetNodesContained().begin();
222 iter!=this->rGetBox(mHalosRight[i]).rGetNodesContained().end();
223 iter++)
224 {
225 mHaloNodesRight.push_back((*iter)->GetIndex());
226 }
227 }
229
230template<unsigned DIM>
232{
233 return mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow();
234}
235
236template<unsigned DIM>
238{
239 return (!(globalIndex<mMinBoxIndex) && !(mMaxBoxIndex<globalIndex));
240}
241
242template<unsigned DIM>
245 bool is_halo_right = ((globalIndex > mMaxBoxIndex) && !(globalIndex > mMaxBoxIndex + mNumBoxesInAFace));
246 bool is_halo_left = ((globalIndex < mMinBoxIndex) && !(globalIndex < mMinBoxIndex - mNumBoxesInAFace));
247
248 // Also need to check for periodic boxes
249 if ( mIsPeriodicAcrossProcs )
250 {
252 {
253 is_halo_right = (globalIndex < mNumBoxesInAFace);
254 }
255 if ( PetscTools::AmMaster() )
257 is_halo_left = (!(globalIndex < (mNumBoxesEachDirection[DIM-1]-1)*mNumBoxesInAFace) && (globalIndex < mNumBoxesEachDirection[DIM-1]*mNumBoxesInAFace ));
258 }
259 }
260
261 return (PetscTools::IsParallel() && (is_halo_right || is_halo_left));
262}
263
264template<unsigned DIM>
267 bool is_on_boundary = !(globalIndex < mMaxBoxIndex - mNumBoxesInAFace) || (globalIndex < mMinBoxIndex + mNumBoxesInAFace);
268
269 return (PetscTools::IsSequential() || !(is_on_boundary));
270}
272template<unsigned DIM>
273unsigned DistributedBoxCollection<DIM>::CalculateGlobalIndex(c_vector<unsigned, DIM> gridIndices)
274{
276 unsigned global_index;
277
278 switch (DIM)
279 {
280 case 1:
282 global_index = gridIndices(0);
283 break;
284 }
285 case 2:
287 global_index = gridIndices(0) +
288 gridIndices(1) * mNumBoxesEachDirection(0);
289 break;
290 }
291 case 3:
292 {
293 global_index = gridIndices(0) +
294 gridIndices(1) * mNumBoxesEachDirection(0) +
295 gridIndices(2) * mNumBoxesEachDirection(0) * mNumBoxesEachDirection(1);
296 break;
297 }
298 default:
299 {
301 }
303 return global_index;
304}
305
306template<unsigned DIM>
308{
309 // Get the location of the node
310 c_vector<double, DIM> location = pNode->rGetLocation();
311 return CalculateContainingBox(location);
312}
313
314template<unsigned DIM>
315unsigned DistributedBoxCollection<DIM>::CalculateContainingBox(c_vector<double, DIM>& rLocation)
317 // The node must lie inside the boundary of the box collection
318 for (unsigned i=0; i<DIM; i++)
319 {
320 if ((rLocation[i] < mDomainSize(2*i)) || !(rLocation[i] < mDomainSize(2*i+1)))
321 {
322 EXCEPTION("The point provided is outside all of the boxes");
324 }
325
326 // Compute the containing box index in each dimension
327 c_vector<unsigned, DIM> containing_box_indices = scalar_vector<unsigned>(DIM, 0u);
328 for (unsigned i=0; i<DIM; i++)
329 {
330 double box_counter = mDomainSize(2*i);
331 while (!((box_counter + mBoxWidth) > rLocation[i] + msFudge))
332 {
333 containing_box_indices[i]++;
334 box_counter += mBoxWidth;
335 }
337
338 // Use these to compute the index of the containing box
339 unsigned containing_box_index = 0;
340 for (unsigned i=0; i<DIM; i++)
341 {
342 unsigned temp = 1;
343 for (unsigned j=0; j<i; j++)
344 {
345 temp *= mNumBoxesEachDirection(j);
346 }
347 containing_box_index += temp*containing_box_indices[i];
349
350 // This index must be less than the total number of boxes
351 assert(containing_box_index < mNumBoxes);
352
353 return containing_box_index;
354}
355
356template<unsigned DIM>
357c_vector<unsigned, DIM> DistributedBoxCollection<DIM>::CalculateGridIndices(unsigned globalIndex)
358{
359 c_vector<unsigned, DIM> grid_indices;
360
361 switch (DIM)
363 case 1:
364 {
365 grid_indices(0) = globalIndex;
366 break;
368 case 2:
369 {
370 unsigned num_x = mNumBoxesEachDirection(0);
371 grid_indices(0) = globalIndex % num_x;
372 grid_indices(1) = (globalIndex - grid_indices(0)) / num_x;
373 break;
375 case 3:
376 {
377 unsigned num_x = mNumBoxesEachDirection(0);
378 unsigned num_xy = mNumBoxesEachDirection(0)*mNumBoxesEachDirection(1);
379 grid_indices(0) = globalIndex % num_x;
380 grid_indices(1) = (globalIndex % num_xy - grid_indices(0)) / num_x;
381 grid_indices(2) = globalIndex / num_xy;
382 break;
383 }
384 default:
387 }
388 }
389
390 return grid_indices;
391}
392
393template<unsigned DIM>
395{
396 // Check first for local ownership
397 if (!(boxIndex<mMinBoxIndex) && !(mMaxBoxIndex<boxIndex))
398 {
399 return mBoxes[boxIndex-mMinBoxIndex];
400 }
402 // If normal execution reaches this point then the box does not belong to the process so we will check for a halo box
403 return rGetHaloBox(boxIndex);
404}
405
406template<unsigned DIM>
408{
409 assert(IsHaloBox(boxIndex));
410
411 unsigned local_index = mHaloBoxesMapping.find(boxIndex)->second;
412
413 return mHaloBoxes[local_index];
414}
415
416template<unsigned DIM>
418{
419 return mNumBoxes;
420}
421
422template<unsigned DIM>
424{
425 return mBoxes.size();
426}
427
428template<unsigned DIM>
430{
431 return mDomainSize;
432}
433
434template<unsigned DIM>
436{
437 return mAreLocalBoxesSet;
438}
439
440template<unsigned DIM>
442{
443 return mBoxWidth;
444}
445
446template<unsigned DIM>
448{
449 return mIsPeriodicInX;
450}
451
452template<unsigned DIM>
454{
455 return mIsPeriodicInY;
456}
457
458template<unsigned DIM>
460{
461 return mIsPeriodicInZ;
462}
463
464template<unsigned DIM>
466{
467 return mIsPeriodicAcrossProcs;
468}
469
470template<unsigned DIM>
472{
473 c_vector<bool, DIM> periodic_dims;
474 periodic_dims(0) = mIsPeriodicInX;
475 if (DIM > 1)
476 {
477 periodic_dims(1) = mIsPeriodicInY;
478 }
479 if (DIM>2)
480 {
481 periodic_dims(2) = mIsPeriodicInZ;
482 }
483 return periodic_dims;
484}
485
486template<unsigned DIM>
488{
489 return mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow();
490}
491
492template<unsigned DIM>
493int DistributedBoxCollection<DIM>::LoadBalance(std::vector<int> localDistribution)
494{
495 MPI_Status status;
496
497 int proc_right = (PetscTools::AmTopMost()) ? MPI_PROC_NULL : (int)PetscTools::GetMyRank() + 1;
498 int proc_left = (PetscTools::AmMaster()) ? MPI_PROC_NULL : (int)PetscTools::GetMyRank() - 1;
499
500 // A variable that will return the new number of rows.
501 int new_rows = localDistribution.size();
502
506 unsigned rows_on_left_process = 0;
507 std::vector<int> node_distr_on_left_process;
508
509 unsigned num_local_rows = localDistribution.size();
510
511 MPI_Send(&num_local_rows, 1, MPI_UNSIGNED, proc_right, 123, PETSC_COMM_WORLD);
512 MPI_Recv(&rows_on_left_process, 1, MPI_UNSIGNED, proc_left, 123, PETSC_COMM_WORLD, &status);
513
514 node_distr_on_left_process.resize(rows_on_left_process > 0 ? rows_on_left_process : 1);
515
516 MPI_Send(&localDistribution[0], num_local_rows, MPI_INT, proc_right, 123, PETSC_COMM_WORLD);
517 MPI_Recv(&node_distr_on_left_process[0], rows_on_left_process, MPI_INT, proc_left, 123, PETSC_COMM_WORLD, &status);
518
522 int local_load = 0;
523 for (unsigned i=0; i<localDistribution.size(); i++)
524 {
525 local_load += localDistribution[i];
526 }
527 int load_on_left_proc = 0;
528 for (unsigned i=0; i<node_distr_on_left_process.size(); i++)
529 {
530 load_on_left_proc += node_distr_on_left_process[i];
531 }
532
534 {
535 // Calculate (Difference in load with a shift) - (Difference in current loads) for a move left and right of the boundary
536 // This code uses integer arithmetic in order to avoid the rounding errors associated with doubles
537 int local_to_left_sq = (local_load - load_on_left_proc) * (local_load - load_on_left_proc);
538 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]) );
539 delta_left = delta_left*delta_left - local_to_left_sq;
540
541 int delta_right = ( (local_load - localDistribution[0]) - (load_on_left_proc + localDistribution[0]));
542 delta_right = delta_right*delta_right - local_to_left_sq;
543
544 // If a delta is negative we should accept that change. If both are negative choose the largest change.
545 int local_change = 0;
546 bool move_left = (!(delta_left > 0) && (node_distr_on_left_process.size() > 1));
547 if (move_left)
548 {
549 local_change = 1;
550 }
551
552 bool move_right = !(delta_right > 0) && (localDistribution.size() > 2);
553 if (move_right)
554 {
555 local_change = -1;
556 }
557
558 if (move_left && move_right)
559 {
560 local_change = (fabs((double)delta_right) > fabs((double)delta_left)) ? -1 : 1;
561 }
562
563 // Update the number of local rows.
564 new_rows += local_change;
565
566 // Send the result of the calculation back to the left processes.
567 MPI_Send(&local_change, 1, MPI_INT, proc_left, 123, PETSC_COMM_WORLD);
568 }
569
570 // Receive changes from right hand process.
571 int remote_change = 0;
572 MPI_Recv(&remote_change, 1, MPI_INT, proc_right, 123, PETSC_COMM_WORLD, &status);
573
574 // Update based on change or right/top boundary
575 new_rows -= remote_change;
576
577 return new_rows;
578}
579
580template<unsigned DIM>
582{
583 if (mAreLocalBoxesSet)
584 {
585 EXCEPTION("Local Boxes Are Already Set");
586 }
587 else
588 {
589 switch (DIM)
590 {
591 case 1:
592 {
593 // We only need to look for neighbours in the current and successive boxes plus some others for halos
594 mLocalBoxes.clear();
595
596 // Iterate over the global box indices
597 for (unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
598 {
599 std::set<unsigned> local_boxes;
600
601 // Insert the current box
602 local_boxes.insert(global_index);
603
604 // Set some bools to find out where we are
605 bool right = (global_index==mNumBoxesEachDirection(0)-1);
606 bool left = (global_index == 0);
607 bool proc_left = (global_index == mpDistributedBoxStackFactory->GetLow());
608
609 // If we're not at the right-most box, then insert the box to the right
610 if (!right)
611 {
612 local_boxes.insert(global_index+1);
613 }
614 // If we're on a left process boundary and not on process 0, insert the (halo) box to the left
615 if (proc_left && !left)
616 {
617 local_boxes.insert(global_index-1);
618 }
619
620 mLocalBoxes.push_back(local_boxes);
621 }
622 break;
623 }
624 case 2:
625 {
626 // Clear the local boxes
627 mLocalBoxes.clear();
628
629 // Now we need to work out the min and max y indices
630 unsigned j_start = CalculateGridIndices(mMinBoxIndex)(1);
631 unsigned j_end = CalculateGridIndices(mMaxBoxIndex)(1);
632 // Determine the number of boxes in each direction
633 unsigned nI = mNumBoxesEachDirection(0);
634 unsigned nJ = mNumBoxesEachDirection(1);
635 // Locally store the bottom processor row
636 unsigned bottom_proc = mpDistributedBoxStackFactory->GetLow();
637 unsigned top_proc = mpDistributedBoxStackFactory->GetHigh();
638
639 // Outer loop: over y
640 for ( unsigned j = j_start; j < (j_end+1); j++ )
641 {
642 // Inner loop: over x
643 for ( unsigned i = 0; i < nI; i++ )
644 {
645 std::set<unsigned> local_boxes;
646 /* We want to add (for non-boundary)
647 (i-1,j+1) (i,j+1) (i+1,j+1)
648 (i,j) (i+1,j) */
649
650 int j_mod = j; // This is used to account for y periodic boundaries
651 // If we are on the bottom of the processor, we may need to add the row below
652 int dj = -1 * (int)(j == bottom_proc && (j > 0 || (mIsPeriodicInY && top_proc < nJ)) );
653 int j_mod_2 = 0;
654 // The min ensures we don't go above the top boundary
655 for (; dj < std::min((int)nJ-j_mod,(int)2); dj++ )
656 {
657 // We need to change to the top row if we are in the condition where dj == -1 and j == 0 which is only
658 // when we are periodic in y and the top row is on a different processor to the bottom row
659 if ( mIsPeriodicInY && j == 0 && dj < 1 && top_proc < nJ )
660 {
661 j_mod_2 = ( dj < 0 ) ? nJ : 0;
662 }
663
664 // The -1*dj ensures we get the upper left, the max ensures we don't hit the left boundary,
665 // the min ensures we don't hit the right boundary
666 int boxi = std::max((int) i-1*std::abs(dj),(int) 0);
667 for ( ; boxi < std::min((int)i+2,(int)nI); boxi++ )
668 {
669 local_boxes.insert( (j_mod+dj+j_mod_2)*nI + boxi );
670 }
671 // Add in the x periodicity at the right boundary, we then want: (0,j) and (0,j+1)
672 if ( i==(nI-1) && mIsPeriodicInX )
673 {
674 local_boxes.insert( (j_mod+dj+j_mod_2)*nI );
675 }
676 // If the y boundary is periodic, the new j level is the bottom row; we use jmod to adjust
677 // for this to adjust for dj = 1
678 if ( mIsPeriodicInY && j == (nJ-1) && dj == 0 )
679 {
680 j_mod = -1;
681 }
682 }
683 // Now add the left-upper box if x periodic and on the boundary
684 if ( mIsPeriodicInX && i == 0 )
685 {
686 if ( j < (nJ-1) )
687 {
688 local_boxes.insert( (j+2)*nI - 1 );
689 if ( j==0 && mIsPeriodicInY && top_proc < nJ )
690 {
691 local_boxes.insert( nI*nJ - 1 );
692 }
693 }
694 else if ( mIsPeriodicInY )
695 {
696 local_boxes.insert( nI-1 );
697 }
698 // If wee are on the bottom of the process need to add
699 // the box on the right in the row below
700 if ( j == bottom_proc && j > 0 )
701 {
702 local_boxes.insert( j*nI - 1 );
703 }
704 }
705
706 // Add to the local boxes
707 mLocalBoxes.push_back(local_boxes);
708 }
709 }
710
711 break;
712 }
713 case 3:
714 {
715 // Clear the local boxes
716 mLocalBoxes.clear();
717
718 // Now we need to work out the min and max z indices
719 unsigned k_start = CalculateGridIndices(mMinBoxIndex)(2);
720 unsigned k_end = CalculateGridIndices(mMaxBoxIndex)(2);
721
722 // Determine the number of boxes in each direction
723 unsigned nI = mNumBoxesEachDirection(0);
724 unsigned nJ = mNumBoxesEachDirection(1);
725 unsigned nK = mNumBoxesEachDirection(2);
726
727 // Work out the bottom/top processor
728 unsigned bottom_proc = mpDistributedBoxStackFactory->GetLow();
729 unsigned top_proc = mpDistributedBoxStackFactory->GetHigh();
730
731 // Outer loop: over z
732 for ( unsigned k = k_start; k <= k_end; k++ )
733 {
734 // Middle loop: over y
735 for ( unsigned j = 0; j < nJ; j++ )
736 {
737 // Inner loop: over x
738 for ( unsigned i = 0; i < nI; i++ )
739 {
740 std::set<unsigned> local_boxes;
741
742 // Same z level
743 unsigned z_offset = k*nI*nJ;
744 // (See case dim=2 for commented code of X-Y implementation)
745 int j_mod = (int)j;
746 for ( int dj = 0; dj < std::min((int)nJ-j_mod,2); dj++ )
747 {
748 for ( int boxi = std::max((int)i-1*dj,0);
749 boxi < std::min((int)i+2,(int)nI); boxi++ )
750 {
751 local_boxes.insert( z_offset + (j_mod+dj)*nI + boxi );
752 }
753 if ( i==(nI-1) && mIsPeriodicInX )
754 {
755 local_boxes.insert( z_offset + (j_mod+dj)*nI );
756 }
757 if ( mIsPeriodicInY && j == (nJ-1) )
758 {
759 j_mod = -1;
760 }
761 }
762 // Now add the left-upper box if x periodic and on the boundary
763 if ( mIsPeriodicInX && i == 0 )
764 {
765 if ( j < (nJ-1) )
766 {
767 local_boxes.insert( z_offset + (j+2)*nI - 1 );
768 }
769 else if ( mIsPeriodicInY )
770 {
771 local_boxes.insert( z_offset + nI-1 );
772 }
773 }
774
775 // Add all the surrounding boxes from the z level above if required
776 std::vector<unsigned> k_offset;
777 if ( k < (nK-1) )
778 {
779 k_offset.push_back(k+1);
780 }
781 if ( k == bottom_proc && k > 0 )
782 {
783 // Need the level below if we are on the lowest processor
784 k_offset.push_back(k-1);
785 }
786 if ( mIsPeriodicInZ && k == (nK-1) )
787 {
788 k_offset.push_back(0);
789 }
790 else if ( mIsPeriodicInZ && k==0 && (top_proc < nK) )
791 {
792 k_offset.push_back(nK-1);
793 }
794
795 for ( std::vector<unsigned>::iterator k_offset_it = k_offset.begin(); k_offset_it != k_offset.end(); ++k_offset_it )
796 {
797 z_offset = (*k_offset_it)*nI*nJ;
798 // Periodicity adjustments
799 int pX = (int) mIsPeriodicInX;
800 int pY = (int) mIsPeriodicInY;
801 for ( int boxi = std::max((int)i-1,-1*pX); boxi < std::min((int)i+2,(int)nI+pX); boxi++ )
802 {
803 for ( int boxj = std::max((int)j-1,-1*pY); boxj < std::min((int)j+2,(int)nJ+pY); boxj++ )
804 {
805 int box_to_add = z_offset + (boxj*(int)nI) + boxi;
806 // Check for x periodicity
807 // i==(nI-1) only when we are periodic and on the right,
808 // i==-1 only when we are periodic and on the left
809 box_to_add += ( (unsigned)(boxi==-1) - (unsigned)(boxi==(int)nI) )*nI;
810 // Check for y periodicity
811 // j==nJ only when we are periodic and on the right,
812 // j==-1 only when we are periodic and on the left
813 box_to_add += ( (unsigned)(boxj==-1) - (unsigned)(boxj==(int)nJ) )*nI*nJ;
814 local_boxes.insert( box_to_add );
815 }
816 }
817 }
818
819 // Add to the local boxes
820 mLocalBoxes.push_back(local_boxes);
821 }
822 }
823 }
824 break;
825 }
826 default:
828 }
829 mAreLocalBoxesSet=true;
830 }
831}
832
833template<unsigned DIM>
835{
836 mAreLocalBoxesSet = true;
837 switch (DIM)
838 {
839 case 1:
840 {
841 for (unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
842 {
843 std::set<unsigned> local_boxes;
844
845 local_boxes.insert(i);
846
847 // add the two neighbours
848 if (i!=0)
849 {
850 local_boxes.insert(i-1);
851 }
852 if (i+1 != mNumBoxesEachDirection(0))
853 {
854 local_boxes.insert(i+1);
855 }
856
857 mLocalBoxes.push_back(local_boxes);
858 }
859 break;
860 }
861 case 2:
862 {
863 mLocalBoxes.clear();
864
865 unsigned M = mNumBoxesEachDirection(0);
866 unsigned N = mNumBoxesEachDirection(1);
867
868 std::vector<bool> is_xmin(N*M); // far left
869 std::vector<bool> is_xmax(N*M); // far right
870 std::vector<bool> is_ymin(N*M); // bottom
871 std::vector<bool> is_ymax(N*M); // top
872
873 for (unsigned i=0; i<M*N; i++)
874 {
875 is_xmin[i] = (i%M==0);
876 is_xmax[i] = ((i+1)%M==0);
877 is_ymin[i] = (i%(M*N)<M);
878 is_ymax[i] = (i%(M*N)>=(N-1)*M);
879 }
880
881 for (unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
882 {
883 std::set<unsigned> local_boxes;
884
885 local_boxes.insert(i);
886
887 // add the box to the left
888 if (!is_xmin[i])
889 {
890 local_boxes.insert(i-1);
891 }
892 else // Add Periodic Box if needed
893 {
894 if (mIsPeriodicInX)
895 {
896 local_boxes.insert(i+M-1);
897 }
898 }
899
900 // add the box to the right
901 if (!is_xmax[i])
902 {
903 local_boxes.insert(i+1);
904 }
905 else // Add Periodic Box if needed
906 {
907 if (mIsPeriodicInX)
908 {
909 local_boxes.insert(i-M+1);
910 }
911 }
912
913 // add the one below
914 if (!is_ymin[i])
915 {
916 local_boxes.insert(i-M);
917 }
918 else // Add periodic box if needed
919 {
920 if (mIsPeriodicInY)
921 {
922 local_boxes.insert(i+(N-1)*M);
923 }
924 }
925
926 // add the one above
927 if (!is_ymax[i])
928 {
929 local_boxes.insert(i+M);
930 }
931 else // Add periodic box if needed
932 {
933 if (mIsPeriodicInY)
934 {
935 local_boxes.insert(i-(N-1)*M);
936 }
937 }
938
939 // add the four corner boxes
940
941 if ( (!is_xmin[i]) && (!is_ymin[i]) )
942 {
943 local_boxes.insert(i-1-M);
944 }
945 if ( (!is_xmin[i]) && (!is_ymax[i]) )
946 {
947 local_boxes.insert(i-1+M);
948 }
949 if ( (!is_xmax[i]) && (!is_ymin[i]) )
950 {
951 local_boxes.insert(i+1-M);
952 }
953 if ( (!is_xmax[i]) && (!is_ymax[i]) )
954 {
955 local_boxes.insert(i+1+M);
956 }
957
958 // Add periodic corner boxes if needed
959 if (mIsPeriodicInX)
960 {
961 if ( (is_xmin[i]) && (!is_ymin[i]) )
962 {
963 local_boxes.insert(i-1);
964 }
965 if ( (is_xmin[i]) && (!is_ymax[i]) )
966 {
967 local_boxes.insert(i-1+2*M);
968 }
969 if ( (is_xmax[i]) && (!is_ymin[i]) )
970 {
971 local_boxes.insert(i+1-2*M);
972 }
973 if ( (is_xmax[i]) && (!is_ymax[i]) )
974 {
975 local_boxes.insert(i+1);
976 }
977 }
978 if(mIsPeriodicInY)
979 {
980 if( (is_ymin[i]) && !(is_xmin[i]) )
981 {
982 local_boxes.insert(i+(N-1)*M-1);
983 }
984 if( (is_ymin[i]) && !(is_xmax[i]) )
985 {
986 local_boxes.insert(i+(N-1)*M+1);
987 }
988 if( (is_ymax[i]) && !(is_xmin[i]) )
989 {
990 local_boxes.insert(i-(N-1)*M-1);
991 }
992 if( (is_ymax[i]) && !(is_xmax[i]) )
993 {
994 local_boxes.insert(i-(N-1)*M+1);
995 }
996 }
997 if(mIsPeriodicInX && mIsPeriodicInY)
998 {
999 if( i==0 ) // Lower left corner
1000 {
1001 local_boxes.insert(M*N-1); // Add upper right corner
1002 }
1003 else if( i==(M-1) ) // Lower right corner
1004 {
1005 local_boxes.insert(M*(N-1)); // Add upper left corner
1006 }
1007 else if( i==(M*(N-1)) ) // Upper left corner
1008 {
1009 local_boxes.insert(M-1); // Add lower right corner
1010 }
1011 else if( i==(M*N-1) ) // Upper right corner
1012 {
1013 local_boxes.insert(0); // Lower left corner
1014 }
1015 }
1016
1017 mLocalBoxes.push_back(local_boxes);
1018 }
1019 break;
1020 }
1021 case 3:
1022 {
1023 mLocalBoxes.clear();
1024
1025 unsigned M = mNumBoxesEachDirection(0);
1026 unsigned N = mNumBoxesEachDirection(1);
1027 unsigned P = mNumBoxesEachDirection(2);
1028
1029 std::vector<bool> is_xmin(N*M*P); // far left
1030 std::vector<bool> is_xmax(N*M*P); // far right
1031 std::vector<bool> is_ymin(N*M*P); // nearest
1032 std::vector<bool> is_ymax(N*M*P); // furthest
1033 std::vector<bool> is_zmin(N*M*P); // bottom layer
1034 std::vector<bool> is_zmax(N*M*P); // top layer
1035
1036 for (unsigned i=0; i<M*N*P; i++)
1037 {
1038 is_xmin[i] = (i%M==0);
1039 is_xmax[i] = ((i+1)%M==0);
1040 is_ymin[i] = (i%(M*N)<M);
1041 is_ymax[i] = (i%(M*N)>=(N-1)*M);
1042 is_zmin[i] = (i<M*N);
1043 is_zmax[i] = (i>=M*N*(P-1));
1044 }
1045
1046 for (unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
1047 {
1048 std::set<unsigned> local_boxes;
1049
1050 // add itself as a local box
1051 local_boxes.insert(i);
1052
1053 // now add all 26 other neighbours.....
1054
1055 // add the box left
1056 if (!is_xmin[i])
1057 {
1058 local_boxes.insert(i-1);
1059
1060 // plus some others towards the left
1061 if (!is_ymin[i])
1062 {
1063 local_boxes.insert(i-1-M);
1064 }
1065
1066 if (!is_ymax[i])
1067 {
1068 local_boxes.insert(i-1+M);
1069 }
1070
1071 if (!is_zmin[i])
1072 {
1073 local_boxes.insert(i-1-M*N);
1074 }
1075
1076 if (!is_zmax[i])
1077 {
1078 local_boxes.insert(i-1+M*N);
1079 }
1080 }
1081
1082 // add the box to the right
1083 if (!is_xmax[i])
1084 {
1085 local_boxes.insert(i+1);
1086
1087 // plus some others towards the right
1088 if (!is_ymin[i])
1089 {
1090 local_boxes.insert(i+1-M);
1091 }
1092
1093 if (!is_ymax[i])
1094 {
1095 local_boxes.insert(i+1+M);
1096 }
1097
1098 if (!is_zmin[i])
1099 {
1100 local_boxes.insert(i+1-M*N);
1101 }
1102
1103 if (!is_zmax[i])
1104 {
1105 local_boxes.insert(i+1+M*N);
1106 }
1107 }
1108
1109 // add the boxes next along the y axis
1110 if (!is_ymin[i])
1111 {
1112 local_boxes.insert(i-M);
1113
1114 // and more in this plane
1115 if (!is_zmin[i])
1116 {
1117 local_boxes.insert(i-M-M*N);
1118 }
1119
1120 if (!is_zmax[i])
1121 {
1122 local_boxes.insert(i-M+M*N);
1123 }
1124 }
1125
1126 // add the boxes next along the y axis
1127 if (!is_ymax[i])
1128 {
1129 local_boxes.insert(i+M);
1130
1131 // and more in this plane
1132 if (!is_zmin[i])
1133 {
1134 local_boxes.insert(i+M-M*N);
1135 }
1136
1137 if (!is_zmax[i])
1138 {
1139 local_boxes.insert(i+M+M*N);
1140 }
1141 }
1142
1143 // add the box directly above
1144 if (!is_zmin[i])
1145 {
1146 local_boxes.insert(i-N*M);
1147 }
1148
1149 // add the box directly below
1150 if (!is_zmax[i])
1151 {
1152 local_boxes.insert(i+N*M);
1153 }
1154
1155 // finally, the 8 corners are left
1156
1157 if ( (!is_xmin[i]) && (!is_ymin[i]) && (!is_zmin[i]) )
1158 {
1159 local_boxes.insert(i-1-M-M*N);
1160 }
1161
1162 if ( (!is_xmin[i]) && (!is_ymin[i]) && (!is_zmax[i]) )
1163 {
1164 local_boxes.insert(i-1-M+M*N);
1165 }
1166
1167 if ( (!is_xmin[i]) && (!is_ymax[i]) && (!is_zmin[i]) )
1168 {
1169 local_boxes.insert(i-1+M-M*N);
1170 }
1171
1172 if ( (!is_xmin[i]) && (!is_ymax[i]) && (!is_zmax[i]) )
1173 {
1174 local_boxes.insert(i-1+M+M*N);
1175 }
1176
1177 if ( (!is_xmax[i]) && (!is_ymin[i]) && (!is_zmin[i]) )
1178 {
1179 local_boxes.insert(i+1-M-M*N);
1180 }
1181
1182 if ( (!is_xmax[i]) && (!is_ymin[i]) && (!is_zmax[i]) )
1183 {
1184 local_boxes.insert(i+1-M+M*N);
1185 }
1186
1187 if ( (!is_xmax[i]) && (!is_ymax[i]) && (!is_zmin[i]) )
1188 {
1189 local_boxes.insert(i+1+M-M*N);
1190 }
1191
1192 if ( (!is_xmax[i]) && (!is_ymax[i]) && (!is_zmax[i]) )
1193 {
1194 local_boxes.insert(i+1+M+M*N);
1195 }
1196
1197 // Now add the periodic boxes if any periodicity
1198 if (mIsPeriodicInX && ( is_xmin[i] || is_xmax[i]) )
1199 {
1200 // We are repeating the same steps on each z level, so easiest is
1201 // to make a vector of the z levels we want
1202 std::vector< int > z_i_offsets(1,0); // Add the current z box level
1203 if ( !is_zmin[i] )
1204 {
1205 z_i_offsets.push_back(-M*N); // Add the z box level below
1206 }
1207 if ( !is_zmax[i] )
1208 {
1209 z_i_offsets.push_back(M*N); // Add the z box level above
1210 }
1211
1212 // If we are on the left, add the nine on the right
1213 if ( is_xmin[i] )
1214 {
1215 // Loop over the z levels
1216 for ( std::vector<int>::iterator it = z_i_offsets.begin(); it != z_i_offsets.end(); it++ )
1217 {
1218 local_boxes.insert( i + (*it) + (M-1) ); // The right-most box on the same row
1219 // We also need to check for y boundaries
1220 if ( !is_ymin[i] )
1221 {
1222 local_boxes.insert( i + (*it) - 1 ); // The right-most box one row below
1223 }
1224 if ( !is_ymax[i] )
1225 {
1226 local_boxes.insert( i + (*it) + (2*M-1) ); // The right-most box one row above
1227 }
1228 }
1229 }
1230
1231 // If we are on the right, add the nine on the left
1232 else if ( is_xmax[i] )
1233 {
1234 // Loop over the z levels
1235 for ( std::vector<int>::iterator it = z_i_offsets.begin(); it != z_i_offsets.end(); it++ )
1236 {
1237 local_boxes.insert( i + (*it) - (M-1) ); // The left-most box on the same row
1238 // We also need to check for y boundaries
1239 if ( !is_ymin[i] )
1240 {
1241 local_boxes.insert( i + (*it) - (2*M-1) ); // The left-most box one row below
1242 }
1243 if ( !is_ymax[i] )
1244 {
1245 local_boxes.insert( i + (*it) + 1 ); // The left-most box one row below
1246 }
1247 }
1248 }
1249 }
1250
1251 if ( mIsPeriodicInY && (is_ymax[i] || is_ymin[i]) )
1252 {
1253 // We consider the current and upper z level and create a vector of the
1254 // opposite box indices that need to be added
1255 std::vector<unsigned> opp_box_i(0);
1256 if ( is_ymin[i] )
1257 {
1258 opp_box_i.push_back(i + (N-1)*M); // Current z level
1259 if ( !is_zmin[i] )
1260 {
1261 opp_box_i.push_back( i - M ); // z level below
1262 }
1263 if ( !is_zmax[i] )
1264 {
1265 opp_box_i.push_back(i + 2*M*N - M); // z level above
1266 }
1267 }
1268 else if ( is_ymax[i] )
1269 {
1270 opp_box_i.push_back( i - (N-1)*M ); // Current z level
1271 if ( !is_zmin[i] )
1272 {
1273 opp_box_i.push_back( i - 2*M*N + M ); // z level below
1274 }
1275 if ( !is_zmax[i] )
1276 {
1277 opp_box_i.push_back( i + M ); // z level above
1278 }
1279 }
1280
1281 // Now we add the different boxes, checking for left and right
1282 for ( std::vector<unsigned>::iterator it_opp_box = opp_box_i.begin(); it_opp_box != opp_box_i.end(); it_opp_box++ )
1283 {
1284 local_boxes.insert( *it_opp_box );
1285 if ( !is_xmin[i] )
1286 {
1287 local_boxes.insert( *it_opp_box - 1 );
1288 }
1289 if ( !is_xmax[i] )
1290 {
1291 local_boxes.insert( *it_opp_box + 1 );
1292 }
1293 }
1294 }
1295
1296 if ( mIsPeriodicInX && mIsPeriodicInY )
1297 {
1298 // Need to add the corners
1299 if ( is_xmin[i] && is_ymin[i] )
1300 {
1301 // Current z level
1302 local_boxes.insert(i+M*N-1);
1303 if ( !is_zmax[i] )
1304 {
1305 // Upper z level
1306 local_boxes.insert(i+2*M*N-1);
1307 }
1308 if ( !is_zmin[i] )
1309 {
1310 // Lower z level
1311 local_boxes.insert(i-1);
1312 }
1313 }
1314 if ( is_xmax[i] && is_ymin[i] )
1315 {
1316 local_boxes.insert(i + (N-2)*M + 1);
1317 if ( !is_zmax[i] )
1318 {
1319 local_boxes.insert(i + M*N + (N-2)*M + 1);
1320 }
1321 if ( !is_zmin[i] )
1322 {
1323 // Lower z level
1324 local_boxes.insert(i-2*M+1);
1325 }
1326 }
1327 if ( is_xmin[i] && is_ymax[i] )
1328 {
1329 local_boxes.insert(i + (N-2)*M - 1);
1330 if (!is_zmax[i])
1331 {
1332 // Upper z level
1333 local_boxes.insert(i - 2*M - 1);
1334 }
1335 if (!is_zmin[i])
1336 {
1337 // Lower z level
1338 local_boxes.insert(i-2*(N-1)*M-1);
1339 }
1340
1341 }
1342 if ( is_xmax[i] && is_ymax[i] )
1343 {
1344 if (!is_zmax[i])
1345 {
1346 // Upper z level
1347 local_boxes.insert(i + 1);
1348 }
1349 if (!is_zmin[i])
1350 {
1351 // Lower z level
1352 local_boxes.insert(i - 2*M*N + 1);
1353 }
1354
1355 }
1356 }
1357
1358 if (mIsPeriodicInZ && (is_zmin[i] || is_zmax[i]))
1359 {
1360 if ( is_zmin[i] )
1361 {
1362 // We need to add the top level
1363 unsigned above_box = i+(P-1)*M*N;
1364 local_boxes.insert(above_box);
1365 if (!is_xmin[i])
1366 {
1367 // Also add the boxes to the left and at the top
1368 local_boxes.insert(above_box-1);
1369 if (!is_ymax[i])
1370 {
1371 local_boxes.insert(above_box+M-1);
1372 }
1373
1374 if (!is_ymin[i])
1375 {
1376 local_boxes.insert(above_box-M-1);
1377 }
1378 }
1379 else if ( mIsPeriodicInX )
1380 {
1381 // Add x periodic box in top layer
1382 local_boxes.insert(above_box+M-1);
1383 if ( !is_ymin[i] )
1384 {
1385 local_boxes.insert(above_box+2*M-1);
1386 }
1387 else if ( mIsPeriodicInY )
1388 {
1389 // Add the xy periodic box in top layer if at y min or max
1390 local_boxes.insert(above_box+M*N-1);
1391 }
1392 if ( !is_ymax[i] )
1393 {
1394 local_boxes.insert(above_box+2*M-1);
1395 }
1396 else if ( mIsPeriodicInY )
1397 {
1398 // Add the xy periodic box in top layer if at y min or max
1399 local_boxes.insert(above_box-M*(N-2)-1);
1400 }
1401 }
1402
1403 if (!is_xmax[i])
1404 {
1405 // Also add the boxes to the left and at the top
1406 local_boxes.insert(above_box+1);
1407 if (!is_ymax[i])
1408 {
1409 local_boxes.insert(above_box+M+1);
1410 }
1411 if (!is_ymin[i])
1412 {
1413 local_boxes.insert(above_box-M+1);
1414 }
1415 }
1416 else if ( mIsPeriodicInX )
1417 {
1418 // Add x periodic box in top layer
1419 local_boxes.insert(above_box - M+1);
1420 if ( mIsPeriodicInY )
1421 {
1422 if ( is_ymin[i] )
1423 {
1424 // Add the xy periodic box in top layer if at y min or max
1425 local_boxes.insert(above_box+M*(N-2)+1);
1426 }
1427 else if ( is_ymax[i] )
1428 {
1429 // Add the xy periodic box in top layer if at y min or max
1430 local_boxes.insert(above_box-M*N+1);
1431 }
1432 }
1433 }
1434
1435
1436 if (!is_ymax[i])
1437 {
1438 local_boxes.insert(above_box+M);
1439 }
1440 else if ( mIsPeriodicInY )
1441 {
1442 local_boxes.insert(above_box-M*(N-1));
1443 if ( !is_xmin[i] )
1444 {
1445 local_boxes.insert(above_box-M*(N-1)-1);
1446 }
1447 if ( !is_xmax[i] )
1448 {
1449 local_boxes.insert(above_box-M*(N-1)+1);
1450 }
1451 }
1452 if (!is_ymin[i])
1453 {
1454 local_boxes.insert(above_box-M);
1455 }
1456 else if ( mIsPeriodicInY )
1457 {
1458 local_boxes.insert(above_box+M*(N-1));
1459 if ( !is_xmin[i] )
1460 {
1461 local_boxes.insert(above_box+M*(N-1)-1);
1462 }
1463 if ( !is_xmax[i] )
1464 {
1465 local_boxes.insert(above_box+M*(N-1)+1);
1466 }
1467 }
1468 }
1469 else if ( is_zmax[i] )
1470 {
1471 // We need to add the bottom level
1472 unsigned below_box = i-(P-1)*M*N;
1473 local_boxes.insert(below_box);
1474 if (!is_xmin[i])
1475 {
1476 // Also add the boxes to the left and at the top
1477 local_boxes.insert(below_box-1);
1478 if (!is_ymax[i])
1479 {
1480 local_boxes.insert(below_box+M-1);
1481 }
1482
1483 if (!is_ymin[i])
1484 {
1485 local_boxes.insert(below_box-M-1);
1486 }
1487 }
1488 else if ( mIsPeriodicInX )
1489 {
1490 // Add x periodic box in top layer
1491 local_boxes.insert(below_box+M-1);
1492 if ( mIsPeriodicInY )
1493 {
1494 if ( is_ymin[i] )
1495 {
1496 // Add the xy periodic box in top layer if at y min or max
1497 local_boxes.insert(below_box+M*N-1);
1498 }
1499 else if ( is_ymax[i] )
1500 {
1501 // Add the xy periodic box in top layer if at y min or max
1502 local_boxes.insert(below_box-M*(N-2)-1);
1503 }
1504 }
1505 }
1506
1507 if (!is_xmax[i])
1508 {
1509 // Also add the boxes to the left and at the top
1510 local_boxes.insert(below_box+1);
1511 if (!is_ymax[i])
1512 {
1513 local_boxes.insert(below_box+M+1);
1514 }
1515 if (!is_ymin[i])
1516 {
1517 local_boxes.insert(below_box-M+1);
1518 }
1519 }
1520 else if ( mIsPeriodicInX )
1521 {
1522 // Add x periodic box in top layer
1523 local_boxes.insert(below_box - M+1);
1524 if ( mIsPeriodicInY )
1525 {
1526 if ( is_ymin[i] )
1527 {
1528 // Add the xy periodic box in top layer if at y min or max
1529 local_boxes.insert(below_box+M*(N-2)+1);
1530 }
1531 else if ( is_ymax[i] )
1532 {
1533 // Add the xy periodic box in top layer if at y min or max
1534 local_boxes.insert(below_box-M*N+1);
1535 }
1536 }
1537 }
1538
1539
1540 if (!is_ymax[i])
1541 {
1542 local_boxes.insert(below_box+M);
1543 }
1544 else if ( mIsPeriodicInY )
1545 {
1546 local_boxes.insert(below_box-M*(N-1));
1547 if ( !is_xmin[i] )
1548 {
1549 local_boxes.insert(below_box-M*(N-1)+1);
1550 }
1551 if (!is_xmax[i])
1552 {
1553 local_boxes.insert(below_box-M*(N-1)-1);
1554 }
1555 }
1556 if (!is_ymin[i])
1557 {
1558 local_boxes.insert(below_box-M);
1559 }
1560 else if ( mIsPeriodicInY )
1561 {
1562 local_boxes.insert(below_box+M*(N-1));
1563 if ( !is_xmin[i] )
1564 {
1565 local_boxes.insert(below_box+M*(N-1)+1);
1566 }
1567 if (!is_xmax[i])
1568 {
1569 local_boxes.insert(below_box+M*(N-1)-1);
1570 }
1571 }
1572 }
1573 }
1574
1575 mLocalBoxes.push_back(local_boxes);
1576 }
1577 break;
1578 }
1579 default:
1581 }
1582}
1583
1584template<unsigned DIM>
1585std::set<unsigned>& DistributedBoxCollection<DIM>::rGetLocalBoxes(unsigned boxIndex)
1586{
1587 // Make sure the box is locally owned
1588 assert(!(boxIndex < mMinBoxIndex) && !(mMaxBoxIndex<boxIndex));
1589 return mLocalBoxes[boxIndex-mMinBoxIndex];
1590}
1591
1592template<unsigned DIM>
1594{
1595 unsigned index = CalculateContainingBox(pNode);
1596
1597 return IsBoxOwned(index);
1598}
1599
1600template<unsigned DIM>
1601bool DistributedBoxCollection<DIM>::IsOwned(c_vector<double, DIM>& location)
1602{
1603 unsigned index = CalculateContainingBox(location);
1604
1605 return IsBoxOwned(index);
1606}
1607
1608template<unsigned DIM>
1610{
1611 unsigned box_index = CalculateContainingBox(pNode);
1612 unsigned containing_process = PetscTools::GetMyRank();
1613
1614 if (box_index > mMaxBoxIndex)
1615 {
1616 containing_process++;
1617 }
1618 else if (box_index < mMinBoxIndex)
1619 {
1620 containing_process--;
1621 }
1622
1623 // Need a special case for periodicity
1624 if ( mIsPeriodicAcrossProcs )
1625 {
1626 if (PetscTools::AmMaster() && box_index > ((mNumBoxesEachDirection[DIM-1]-1)*mNumBoxesInAFace-1))
1627 {
1628 // It needs to move to the top process
1629 containing_process = PetscTools::GetNumProcs()-1;
1630 }
1631 else if (PetscTools::AmTopMost() && box_index < mNumBoxesInAFace)
1632 {
1633 // It needs to move to the bottom process
1634 containing_process = 0;
1635 }
1636 }
1637
1638 return containing_process;
1639}
1640
1641template<unsigned DIM>
1643{
1644 return mHaloNodesRight;
1645}
1646
1647template<unsigned DIM>
1649{
1650 return mHaloNodesLeft;
1651}
1652
1653template<unsigned DIM>
1655{
1656 mCalculateNodeNeighbours = calculateNodeNeighbours;
1657}
1658
1659template<unsigned DIM>
1660void DistributedBoxCollection<DIM>::CalculateNodePairs(std::vector<Node<DIM>*>& rNodes, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs)
1661{
1662 rNodePairs.clear();
1663
1664 // Create an empty neighbours set for each node
1665 for (unsigned i=0; i<rNodes.size(); i++)
1666 {
1667 // Get the box containing this node as only nodes on this process have NodeAttributes
1668 // and therefore Neighbours setup.
1669 unsigned box_index = CalculateContainingBox(rNodes[i]);
1670
1671 if (IsBoxOwned(box_index))
1672 {
1673 rNodes[i]->ClearNeighbours();
1674 }
1675 }
1676
1677 for (unsigned box_index=mMinBoxIndex; box_index<=mMaxBoxIndex; box_index++)
1678 {
1679 AddPairsFromBox(box_index, rNodePairs);
1680 }
1681
1682 if (mCalculateNodeNeighbours)
1683 {
1684 for (unsigned i = 0; i < rNodes.size(); i++)
1685 {
1686 // Get the box containing this node as only nodes on this process have NodeAttributes
1687 // and therefore Neighbours setup.
1688 unsigned box_index = CalculateContainingBox(rNodes[i]);
1689
1690 if (IsBoxOwned(box_index))
1691 {
1692 rNodes[i]->RemoveDuplicateNeighbours();
1693 }
1694 }
1695 }
1696}
1697
1698template<unsigned DIM>
1699void DistributedBoxCollection<DIM>::CalculateInteriorNodePairs(std::vector<Node<DIM>*>& rNodes, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs)
1700{
1701 rNodePairs.clear();
1702
1703 // Create an empty neighbours set for each node
1704 for (unsigned i=0; i<rNodes.size(); i++)
1705 {
1706 // Get the box containing this node as only nodes on this process have NodeAttributes
1707 // and therefore Neighbours setup.
1708 unsigned box_index = CalculateContainingBox(rNodes[i]);
1709
1710 if (IsBoxOwned(box_index))
1711 {
1712 rNodes[i]->ClearNeighbours();
1713 rNodes[i]->SetNeighboursSetUp(false);
1714 }
1715 }
1716
1717 for (unsigned box_index=mMinBoxIndex; box_index<=mMaxBoxIndex; box_index++)
1718 {
1719 if (IsInteriorBox(box_index))
1720 {
1721 AddPairsFromBox(box_index, rNodePairs);
1722 }
1723 }
1724
1725 if (mCalculateNodeNeighbours)
1726 {
1727 for (unsigned i = 0; i < rNodes.size(); i++)
1728 {
1729 // Get the box containing this node as only nodes on this process have NodeAttributes
1730 // and therefore Neighbours setup.
1731 unsigned box_index = CalculateContainingBox(rNodes[i]);
1732
1733 if (IsBoxOwned(box_index))
1734 {
1735 rNodes[i]->RemoveDuplicateNeighbours();
1736 rNodes[i]->SetNeighboursSetUp(true);
1737 }
1738 }
1739 }
1740}
1741
1742template<unsigned DIM>
1743void DistributedBoxCollection<DIM>::CalculateBoundaryNodePairs(std::vector<Node<DIM>*>& rNodes, std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs)
1744{
1745 for (unsigned box_index=mMinBoxIndex; box_index<=mMaxBoxIndex; box_index++)
1746 {
1747 if (!IsInteriorBox(box_index))
1748 {
1749 AddPairsFromBox(box_index, rNodePairs);
1750 }
1751 }
1752
1753 if (mCalculateNodeNeighbours)
1754 {
1755 for (unsigned i = 0; i < rNodes.size(); i++)
1756 {
1757 // Get the box containing this node as only nodes on this process have NodeAttributes
1758 // and therefore Neighbours setup.
1759 unsigned box_index = CalculateContainingBox(rNodes[i]);
1760
1761 if (IsBoxOwned(box_index))
1762 {
1763 rNodes[i]->RemoveDuplicateNeighbours();
1764 rNodes[i]->SetNeighboursSetUp(true);
1765 }
1766 }
1767 }
1768}
1769
1770template<unsigned DIM>
1772 std::vector<std::pair<Node<DIM>*, Node<DIM>*> >& rNodePairs)
1773{
1774 // Get the box
1775 Box<DIM>& r_box = rGetBox(boxIndex);
1776
1777 // Get the set of nodes in this box
1778 const std::set< Node<DIM>* >& r_contained_nodes = r_box.rGetNodesContained();
1779
1780 // Get the local boxes to this box
1781 const std::set<unsigned>& local_boxes_indices = rGetLocalBoxes(boxIndex);
1782
1783 // Loop over all the local boxes
1784 for (std::set<unsigned>::iterator box_iter = local_boxes_indices.begin();
1785 box_iter != local_boxes_indices.end();
1786 box_iter++)
1787 {
1788 Box<DIM>* p_neighbour_box;
1789
1790 // Establish whether box is locally owned or halo.
1791 if (IsBoxOwned(*box_iter))
1792 {
1793 p_neighbour_box = &mBoxes[*box_iter - mMinBoxIndex];
1794 }
1795 else // Assume it is a halo.
1796 {
1797 p_neighbour_box = &mHaloBoxes[mHaloBoxesMapping[*box_iter]];
1798 }
1799 assert(p_neighbour_box);
1800
1801 // Get the set of nodes contained in this box
1802 std::set< Node<DIM>* >& r_contained_neighbour_nodes = p_neighbour_box->rGetNodesContained();
1803
1804 // Loop over these nodes
1805 for (typename std::set<Node<DIM>*>::iterator neighbour_node_iter = r_contained_neighbour_nodes.begin();
1806 neighbour_node_iter != r_contained_neighbour_nodes.end();
1807 ++neighbour_node_iter)
1808 {
1809 // Get the index of the other node
1810 unsigned other_node_index = (*neighbour_node_iter)->GetIndex();
1811
1812 // Loop over nodes in this box
1813 for (typename std::set<Node<DIM>*>::iterator node_iter = r_contained_nodes.begin();
1814 node_iter != r_contained_nodes.end();
1815 ++node_iter)
1816 {
1817 unsigned node_index = (*node_iter)->GetIndex();
1818
1819 // If we're in the same box, then take care not to store the node pair twice
1820 if (*box_iter != boxIndex || other_node_index > node_index)
1821 {
1822 rNodePairs.push_back(std::pair<Node<DIM>*, Node<DIM>*>((*node_iter), (*neighbour_node_iter)));
1823 if (mCalculateNodeNeighbours)
1824 {
1825 (*node_iter)->AddNeighbour(other_node_index);
1826 (*neighbour_node_iter)->AddNeighbour(node_index);
1827 }
1828 }
1829
1830 }
1831 }
1832 }
1833}
1834
1835template<unsigned DIM>
1837{
1838 std::vector<int> cell_numbers(mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow(), 0);
1839
1840 for (unsigned global_index=mMinBoxIndex; global_index<=mMaxBoxIndex; global_index++)
1841 {
1842 c_vector<unsigned, DIM> coords = CalculateGridIndices(global_index);
1843 unsigned location_in_vector = coords[DIM-1] - mpDistributedBoxStackFactory->GetLow();
1844 unsigned local_index = global_index - mMinBoxIndex;
1845 cell_numbers[location_in_vector] += mBoxes[local_index].rGetNodesContained().size();
1846 }
1847
1848 return cell_numbers;
1849}
1850
1852
1853template class DistributedBoxCollection<1>;
1854template class DistributedBoxCollection<2>;
1855template class DistributedBoxCollection<3>;
#define EXCEPTION(message)
#define NEVER_REACHED
Definition Box.hpp:51
std::set< Node< DIM > * > & rGetNodesContained()
Definition Box.cpp:56
void AddPairsFromBox(unsigned boxIndex, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs)
std::vector< Box< DIM > > mBoxes
std::set< unsigned > & rGetLocalBoxes(unsigned boxIndex)
DistributedBoxCollection(double boxWidth, c_vector< double, 2 *DIM > domainSize, bool isPeriodicInX=false, bool isPeriodicInY=false, bool isPeriodicInZ=false, int localRows=PETSC_DECIDE)
Box< DIM > & rGetHaloBox(unsigned boxIndex)
c_vector< double, 2 *DIM > rGetDomainSize() const
int LoadBalance(std::vector< int > localDistribution)
void CalculateNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs)
void CalculateInteriorNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs)
unsigned GetProcessOwningNode(Node< DIM > *pNode)
void SetCalculateNodeNeighbours(bool calculateNodeNeighbours)
bool IsInteriorBox(unsigned globalIndex)
unsigned CalculateContainingBox(Node< DIM > *pNode)
std::vector< unsigned > & rGetHaloNodesLeft()
std::vector< int > CalculateNumberOfNodesInEachStrip()
bool IsBoxOwned(unsigned globalIndex)
std::vector< unsigned > & rGetHaloNodesRight()
bool IsHaloBox(unsigned globalIndex)
DistributedVectorFactory * mpDistributedBoxStackFactory
c_vector< bool, DIM > GetIsPeriodicAllDims() const
c_vector< double, 2 *DIM > mDomainSize
void CalculateBoundaryNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs)
c_vector< unsigned, DIM > mNumBoxesEachDirection
unsigned CalculateGlobalIndex(c_vector< unsigned, DIM > gridIndices)
c_vector< unsigned, DIM > CalculateGridIndices(unsigned globalIndex)
Box< DIM > & rGetBox(unsigned boxIndex)
Definition Node.hpp:59
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition Node.cpp:139
static bool AmMaster()
static bool AmTopMost()
static bool IsParallel()
static bool IsSequential()
static unsigned GetMyRank()
static unsigned GetNumProcs()