35 #include "DistributedBoxCollection.hpp"
38 #include "Warnings.hpp"
41 template<
unsigned DIM>
44 template<
unsigned DIM>
46 : mBoxWidth(boxWidth),
47 mIsPeriodicInX(isPeriodicInX),
48 mAreLocalBoxesSet(false),
49 mCalculateNodeNeighbours(true)
58 for (
unsigned i=0; i<DIM; i++)
60 double r = fmod((domainSize[2*i+1]-domainSize[2*i]), boxWidth);
63 domainSize[2*i+1] += boxWidth - r;
72 for (
unsigned i=0; i<DIM; i++)
85 WARNING(
"There are more processes than convenient for the domain/mesh/box size. The domain size has been swollen.")
95 for (
unsigned i=1; i<DIM; i++)
110 c_vector<double, 2*DIM> arbitrary_location;
111 for (
unsigned i=0; i<num_boxes; i++)
113 Box<DIM> new_box(arbitrary_location);
114 mBoxes.push_back(new_box);
116 unsigned global_index = mMinBoxIndex + i;
123 template<
unsigned DIM>
126 delete mpDistributedBoxStackFactory;
129 template<
unsigned DIM>
132 for (
unsigned i=0; i<mBoxes.size(); i++)
134 mBoxes[i].ClearNodes();
136 for (
unsigned i=0; i<mHaloBoxes.size(); i++)
138 mHaloBoxes[i].ClearNodes();
142 template<
unsigned DIM>
146 unsigned Hi = mpDistributedBoxStackFactory->GetHigh();
147 unsigned Lo = mpDistributedBoxStackFactory->GetLow();
152 for (
unsigned i=0; i < mNumBoxesInAFace; i++)
154 c_vector<double, 2*DIM> arbitrary_location;
155 Box<DIM> new_box(arbitrary_location);
156 mHaloBoxes.push_back(new_box);
158 unsigned global_index = Hi * mNumBoxesInAFace + i;
159 mHaloBoxesMapping[global_index] = mHaloBoxes.size()-1;
160 mHalosRight.push_back(global_index - mNumBoxesInAFace);
167 for (
unsigned i=0; i< mNumBoxesInAFace; i++)
169 c_vector<double, 2*DIM> arbitrary_location;
170 Box<DIM> new_box(arbitrary_location);
171 mHaloBoxes.push_back(new_box);
173 unsigned global_index = (Lo - 1) * mNumBoxesInAFace + i;
174 mHaloBoxesMapping[global_index] = mHaloBoxes.size() - 1;
176 mHalosLeft.push_back(global_index + mNumBoxesInAFace);
181 template<
unsigned DIM>
184 mHaloNodesLeft.clear();
185 for (
unsigned i=0; i<mHalosLeft.size(); i++)
187 for (
typename std::set<
Node<DIM>* >::iterator iter=this->rGetBox(mHalosLeft[i]).rGetNodesContained().begin();
188 iter!=this->rGetBox(mHalosLeft[i]).rGetNodesContained().end();
191 mHaloNodesLeft.push_back((*iter)->GetIndex());
196 mHaloNodesRight.clear();
197 for (
unsigned i=0; i<mHalosRight.size(); i++)
199 for (
typename std::set<
Node<DIM>* >::iterator iter=this->rGetBox(mHalosRight[i]).rGetNodesContained().begin();
200 iter!=this->rGetBox(mHalosRight[i]).rGetNodesContained().end();
203 mHaloNodesRight.push_back((*iter)->GetIndex());
208 template<
unsigned DIM>
211 return mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow();
214 template<
unsigned DIM>
217 return (!(globalIndex<mMinBoxIndex) && !(mMaxBoxIndex<globalIndex));
220 template<
unsigned DIM>
223 bool is_halo_right = ((globalIndex > mMaxBoxIndex) && !(globalIndex > mMaxBoxIndex + mNumBoxesInAFace));
224 bool is_halo_left = ((globalIndex < mMinBoxIndex) && !(globalIndex < mMinBoxIndex - mNumBoxesInAFace));
229 template<
unsigned DIM>
232 bool is_on_boundary = !(globalIndex < mMaxBoxIndex - mNumBoxesInAFace) || (globalIndex < mMinBoxIndex + mNumBoxesInAFace);
237 template<
unsigned DIM>
240 unsigned containing_box_index = 0;
241 for (
unsigned i=0; i<DIM; i++)
244 for (
unsigned j=0; j<i; j++)
246 temp *= mNumBoxesEachDirection(j);
248 containing_box_index += temp*coordinateIndices[i];
251 return containing_box_index;
254 template<
unsigned DIM>
259 return CalculateContainingBox(location);
263 template<
unsigned DIM>
267 for (
unsigned i=0; i<DIM; i++)
269 if ( (rLocation[i] < mDomainSize(2*i)) || !(rLocation[i] < mDomainSize(2*i+1)) )
271 EXCEPTION(
"The point provided is outside all of the boxes");
276 c_vector<unsigned, DIM> containing_box_indices = scalar_vector<unsigned>(DIM, 0u);
277 for (
unsigned i=0; i<DIM; i++)
279 double box_counter = mDomainSize(2*i);
280 while (!((box_counter + mBoxWidth) > rLocation[i] + msFudge))
282 containing_box_indices[i]++;
283 box_counter += mBoxWidth;
288 unsigned containing_box_index = 0;
289 for (
unsigned i=0; i<DIM; i++)
292 for (
unsigned j=0; j<i; j++)
294 temp *= mNumBoxesEachDirection(j);
296 containing_box_index += temp*containing_box_indices[i];
300 assert(containing_box_index < mNumBoxes);
302 return containing_box_index;
305 template<
unsigned DIM>
308 c_vector<unsigned, DIM> indices;
314 indices[0]=globalIndex;
319 unsigned remainder=globalIndex % mNumBoxesEachDirection(0);
320 indices[0]=remainder;
321 indices[1]=(
unsigned)(globalIndex/mNumBoxesEachDirection(0));
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)));
343 template<
unsigned DIM>
346 assert(!(boxIndex<mMinBoxIndex) && !(mMaxBoxIndex<boxIndex));
347 return mBoxes[boxIndex-mMinBoxIndex];
350 template<
unsigned DIM>
353 assert(GetHaloBoxOwnership(boxIndex));
355 unsigned local_index = mHaloBoxesMapping.find(boxIndex)->second;
357 return mHaloBoxes[local_index];
360 template<
unsigned DIM>
366 template<
unsigned DIM>
369 return mBoxes.size();
372 template<
unsigned DIM>
378 template<
unsigned DIM>
381 return mAreLocalBoxesSet;
384 template<
unsigned DIM>
390 template<
unsigned DIM>
393 return mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow();
396 template<
unsigned DIM>
405 int new_rows = localDistribution.size();
410 unsigned rows_on_left_process = 0;
411 std::vector<int> node_distr_on_left_process;
413 unsigned num_local_rows = localDistribution.size();
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);
418 node_distr_on_left_process.resize(rows_on_left_process > 0 ? rows_on_left_process : 1);
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);
427 for (
unsigned i=0; i<localDistribution.size(); i++)
429 local_load += localDistribution[i];
431 int load_on_left_proc = 0;
432 for (
unsigned i=0; i<node_distr_on_left_process.size(); i++)
434 load_on_left_proc += node_distr_on_left_process[i];
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;
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;
449 int local_change = 0;
450 bool move_left = (!(delta_left > 0) && (node_distr_on_left_process.size() > 1));
456 bool move_right = !(delta_right > 0) && (localDistribution.size() > 2);
462 if (move_left && move_right)
464 local_change = (fabs((
double)delta_right) > fabs((
double)delta_left)) ? -1 : 1;
468 new_rows += local_change;
471 MPI_Send(&local_change, 1, MPI_INT, proc_left, 123, PETSC_COMM_WORLD);
475 int remote_change = 0;
476 MPI_Recv(&remote_change, 1, MPI_INT, proc_right, 123, PETSC_COMM_WORLD, &status);
479 new_rows -= remote_change;
484 template<
unsigned DIM>
487 if (mAreLocalBoxesSet)
489 EXCEPTION(
"Local Boxes Are Already Set");
501 for (
unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
503 std::set<unsigned> local_boxes;
506 local_boxes.insert(global_index);
509 bool right = (global_index==mNumBoxesEachDirection(0)-1);
510 bool left = (global_index == 0);
511 bool proc_left = (global_index == mpDistributedBoxStackFactory->GetLow());
516 local_boxes.insert(global_index+1);
519 if (proc_left && !left)
521 local_boxes.insert(global_index-1);
524 mLocalBoxes.push_back(local_boxes);
533 for (
unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
535 std::set<unsigned> local_boxes;
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());
545 local_boxes.insert(global_index);
548 if(!bottom && bottom_proc)
550 local_boxes.insert(global_index - mNumBoxesEachDirection(0));
553 local_boxes.insert(global_index - mNumBoxesEachDirection(0) - 1);
557 local_boxes.insert(global_index - mNumBoxesEachDirection(0) + 1);
564 local_boxes.insert(global_index + mNumBoxesEachDirection(0));
568 local_boxes.insert(global_index + mNumBoxesEachDirection(0) + 1);
572 local_boxes.insert(global_index + mNumBoxesEachDirection(0) - 1);
574 else if ( (global_index % mNumBoxesEachDirection(0) == 0) && (mIsPeriodicInX) )
576 local_boxes.insert(global_index + 2 * mNumBoxesEachDirection(0) - 1);
583 local_boxes.insert(global_index + 1);
586 else if ( (global_index % mNumBoxesEachDirection(0) == mNumBoxesEachDirection(0)-1) && (mIsPeriodicInX) )
588 local_boxes.insert(global_index - mNumBoxesEachDirection(0) + 1);
590 if (global_index < mBoxes.size() - mNumBoxesEachDirection(0))
592 local_boxes.insert(global_index + 1);
596 mLocalBoxes.push_back(local_boxes);
604 unsigned num_boxes_xy = mNumBoxesEachDirection(0)*mNumBoxesEachDirection(1);
607 for (
unsigned global_index = mMinBoxIndex; global_index<mMaxBoxIndex+1; global_index++)
609 std::set<unsigned> local_boxes;
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);
622 local_boxes.insert(global_index);
630 local_boxes.insert( global_index - num_boxes_xy + mNumBoxesEachDirection(0) );
633 local_boxes.insert( global_index - num_boxes_xy + mNumBoxesEachDirection(0) - 1);
637 local_boxes.insert( global_index - num_boxes_xy + mNumBoxesEachDirection(0) + 1);
642 local_boxes.insert( global_index - num_boxes_xy + 1);
648 local_boxes.insert( global_index - num_boxes_xy );
652 local_boxes.insert( global_index - num_boxes_xy - 1);
656 local_boxes.insert( global_index - num_boxes_xy - mNumBoxesEachDirection(0));
660 local_boxes.insert( global_index - num_boxes_xy - mNumBoxesEachDirection(0) - 1);
664 local_boxes.insert( global_index - num_boxes_xy - mNumBoxesEachDirection(0) + 1);
672 local_boxes.insert( global_index + 1);
677 local_boxes.insert( global_index + mNumBoxesEachDirection(0));
681 local_boxes.insert( global_index + mNumBoxesEachDirection(0) + 1);
685 local_boxes.insert( global_index + mNumBoxesEachDirection(0) - 1);
692 local_boxes.insert(global_index + num_boxes_xy);
696 local_boxes.insert(global_index + num_boxes_xy + 1);
700 local_boxes.insert(global_index + num_boxes_xy + mNumBoxesEachDirection(0));
703 local_boxes.insert(global_index + num_boxes_xy + mNumBoxesEachDirection(0) + 1);
707 local_boxes.insert(global_index + num_boxes_xy + mNumBoxesEachDirection(0) - 1);
715 local_boxes.insert(global_index + num_boxes_xy - 1);
719 local_boxes.insert(global_index + num_boxes_xy - mNumBoxesEachDirection(0));
723 local_boxes.insert(global_index + num_boxes_xy - mNumBoxesEachDirection(0) - 1);
727 local_boxes.insert(global_index + num_boxes_xy - mNumBoxesEachDirection(0) + 1);
732 mLocalBoxes.push_back(local_boxes);
740 mAreLocalBoxesSet=
true;
746 template<
unsigned DIM>
749 mAreLocalBoxesSet =
true;
754 for (
unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
756 std::set<unsigned> local_boxes;
758 local_boxes.insert(i);
763 local_boxes.insert(i-1);
765 if (i+1 != mNumBoxesEachDirection(0))
767 local_boxes.insert(i+1);
770 mLocalBoxes.push_back(local_boxes);
778 unsigned M = mNumBoxesEachDirection(0);
779 unsigned N = mNumBoxesEachDirection(1);
781 std::vector<bool> is_xmin(N*M);
782 std::vector<bool> is_xmax(N*M);
783 std::vector<bool> is_ymin(N*M);
784 std::vector<bool> is_ymax(N*M);
786 for (
unsigned i=0; i<M*N; i++)
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);
794 for (
unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
796 std::set<unsigned> local_boxes;
798 local_boxes.insert(i);
803 local_boxes.insert(i-1);
809 local_boxes.insert(i+M-1);
816 local_boxes.insert(i+1);
822 local_boxes.insert(i-M+1);
829 local_boxes.insert(i-M);
835 local_boxes.insert(i+M);
840 if ( (!is_xmin[i]) && (!is_ymin[i]) )
842 local_boxes.insert(i-1-M);
844 if ( (!is_xmin[i]) && (!is_ymax[i]) )
846 local_boxes.insert(i-1+M);
848 if ( (!is_xmax[i]) && (!is_ymin[i]) )
850 local_boxes.insert(i+1-M);
852 if ( (!is_xmax[i]) && (!is_ymax[i]) )
854 local_boxes.insert(i+1+M);
860 if( (is_xmin[i]) && (!is_ymin[i]) )
862 local_boxes.insert(i-1);
864 if ( (is_xmin[i]) && (!is_ymax[i]) )
866 local_boxes.insert(i-1+2*M);
868 if ( (is_xmax[i]) && (!is_ymin[i]) )
870 local_boxes.insert(i+1-2*M);
872 if ( (is_xmax[i]) && (!is_ymax[i]) )
874 local_boxes.insert(i+1);
878 mLocalBoxes.push_back(local_boxes);
886 unsigned M = mNumBoxesEachDirection(0);
887 unsigned N = mNumBoxesEachDirection(1);
888 unsigned P = mNumBoxesEachDirection(2);
890 std::vector<bool> is_xmin(N*M*P);
891 std::vector<bool> is_xmax(N*M*P);
892 std::vector<bool> is_ymin(N*M*P);
893 std::vector<bool> is_ymax(N*M*P);
894 std::vector<bool> is_zmin(N*M*P);
895 std::vector<bool> is_zmax(N*M*P);
897 for (
unsigned i=0; i<M*N*P; i++)
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));
907 for (
unsigned i=mMinBoxIndex; i<mMaxBoxIndex+1; i++)
909 std::set<unsigned> local_boxes;
912 local_boxes.insert(i);
919 local_boxes.insert(i-1);
924 local_boxes.insert(i-1-M);
929 local_boxes.insert(i-1+M);
934 local_boxes.insert(i-1-M*N);
939 local_boxes.insert(i-1+M*N);
946 local_boxes.insert(i+1);
951 local_boxes.insert(i+1-M);
956 local_boxes.insert(i+1+M);
961 local_boxes.insert(i+1-M*N);
966 local_boxes.insert(i+1+M*N);
973 local_boxes.insert(i-M);
978 local_boxes.insert(i-M-M*N);
983 local_boxes.insert(i-M+M*N);
990 local_boxes.insert(i+M);
995 local_boxes.insert(i+M-M*N);
1000 local_boxes.insert(i+M+M*N);
1007 local_boxes.insert(i-N*M);
1013 local_boxes.insert(i+N*M);
1018 if ( (!is_xmin[i]) && (!is_ymin[i]) && (!is_zmin[i]) )
1020 local_boxes.insert(i-1-M-M*N);
1023 if ( (!is_xmin[i]) && (!is_ymin[i]) && (!is_zmax[i]) )
1025 local_boxes.insert(i-1-M+M*N);
1028 if ( (!is_xmin[i]) && (!is_ymax[i]) && (!is_zmin[i]) )
1030 local_boxes.insert(i-1+M-M*N);
1033 if ( (!is_xmin[i]) && (!is_ymax[i]) && (!is_zmax[i]) )
1035 local_boxes.insert(i-1+M+M*N);
1038 if ( (!is_xmax[i]) && (!is_ymin[i]) && (!is_zmin[i]) )
1040 local_boxes.insert(i+1-M-M*N);
1043 if ( (!is_xmax[i]) && (!is_ymin[i]) && (!is_zmax[i]) )
1045 local_boxes.insert(i+1-M+M*N);
1048 if ( (!is_xmax[i]) && (!is_ymax[i]) && (!is_zmin[i]) )
1050 local_boxes.insert(i+1+M-M*N);
1053 if ( (!is_xmax[i]) && (!is_ymax[i]) && (!is_zmax[i]) )
1055 local_boxes.insert(i+1+M+M*N);
1058 mLocalBoxes.push_back(local_boxes);
1067 template<
unsigned DIM>
1071 assert(!(boxIndex < mMinBoxIndex) && !(mMaxBoxIndex<boxIndex));
1072 return mLocalBoxes[boxIndex-mMinBoxIndex];
1075 template<
unsigned DIM>
1078 unsigned index = CalculateContainingBox(pNode);
1080 return GetBoxOwnership(index);
1083 template<
unsigned DIM>
1086 unsigned index = CalculateContainingBox(location);
1088 return GetBoxOwnership(index);
1091 template<
unsigned DIM>
1094 unsigned box_index = CalculateContainingBox(pNode);
1097 if (box_index > mMaxBoxIndex)
1099 containing_process++;
1101 else if (box_index < mMinBoxIndex)
1103 containing_process--;
1106 return containing_process;
1109 template<
unsigned DIM>
1112 return mHaloNodesRight;
1115 template<
unsigned DIM>
1118 return mHaloNodesLeft;
1121 template<
unsigned DIM>
1124 mCalculateNodeNeighbours = calculateNodeNeighbours;
1127 template<
unsigned DIM>
1131 rNodeNeighbours.clear();
1134 for (
unsigned i=0; i<rNodes.size(); i++)
1136 unsigned node_index = rNodes[i]->GetIndex();
1139 unsigned box_index = CalculateContainingBox(rNodes[i]);
1141 if (GetBoxOwnership(box_index))
1143 rNodeNeighbours[node_index] = std::set<unsigned>();
1147 for (std::map<unsigned, unsigned>::iterator map_iter = mBoxesMapping.begin();
1148 map_iter != mBoxesMapping.end();
1152 unsigned box_index = map_iter->first;
1154 AddPairsFromBox(box_index, rNodePairs, rNodeNeighbours);
1158 template<
unsigned DIM>
1162 rNodeNeighbours.clear();
1165 for (
unsigned i=0; i<rNodes.size(); i++)
1167 unsigned node_index = rNodes[i]->GetIndex();
1170 unsigned box_index = CalculateContainingBox(rNodes[i]);
1172 if (GetBoxOwnership(box_index))
1174 rNodeNeighbours[node_index] = std::set<unsigned>();
1178 for (std::map<unsigned, unsigned>::iterator map_iter = mBoxesMapping.begin();
1179 map_iter != mBoxesMapping.end();
1183 unsigned box_index = map_iter->first;
1185 if (IsInteriorBox(box_index))
1187 AddPairsFromBox(box_index, rNodePairs, rNodeNeighbours);
1192 template<
unsigned DIM>
1195 for (std::map<unsigned, unsigned>::iterator map_iter = mBoxesMapping.begin();
1196 map_iter != mBoxesMapping.end();
1200 unsigned box_index = map_iter->first;
1202 if (!IsInteriorBox(box_index))
1204 AddPairsFromBox(box_index, rNodePairs, rNodeNeighbours);
1209 template<
unsigned DIM>
1219 std::set<unsigned> local_boxes_indices = GetLocalBoxes(boxIndex);
1222 for (std::set<unsigned>::iterator box_iter = local_boxes_indices.begin();
1223 box_iter != local_boxes_indices.end();
1229 if (GetBoxOwnership(*box_iter))
1231 p_neighbour_box = &mBoxes[mBoxesMapping[*box_iter]];
1235 p_neighbour_box = &mHaloBoxes[mHaloBoxesMapping[*box_iter]];
1237 assert(p_neighbour_box);
1240 std::set< Node<DIM>* >& r_contained_neighbour_nodes = p_neighbour_box->
rGetNodesContained();
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)
1248 unsigned other_node_index = (*neighbour_node_iter)->GetIndex();
1251 for (
typename std::set<
Node<DIM>*>::iterator node_iter = r_contained_nodes.begin();
1252 node_iter != r_contained_nodes.end();
1255 unsigned node_index = (*node_iter)->GetIndex();
1258 if (*box_iter == boxIndex)
1260 if (other_node_index > node_index)
1262 rNodePairs.push_back(std::pair<
Node<DIM>*,
Node<DIM>*>((*node_iter), (*neighbour_node_iter)));
1263 if (mCalculateNodeNeighbours)
1265 rNodeNeighbours[node_index].insert(other_node_index);
1266 rNodeNeighbours[other_node_index].insert(node_index);
1272 rNodePairs.push_back(std::pair<
Node<DIM>*,
Node<DIM>*>((*node_iter), (*neighbour_node_iter)));
1273 if (mCalculateNodeNeighbours)
1275 rNodeNeighbours[node_index].insert(other_node_index);
1276 rNodeNeighbours[other_node_index].insert(node_index);
1285 template<
unsigned DIM>
1288 std::vector<int> cell_numbers(mpDistributedBoxStackFactory->GetHigh() - mpDistributedBoxStackFactory->GetLow(), 0);
1290 for (std::map<unsigned, unsigned>::iterator iter = mBoxesMapping.begin();
1291 iter != mBoxesMapping.end();
1294 c_vector<unsigned, DIM> coords = CalculateCoordinateIndices(iter->first);
1295 unsigned location_in_vector = coords[DIM-1]-mpDistributedBoxStackFactory->GetLow();
1297 cell_numbers[location_in_vector] += mBoxes[iter->second].rGetNodesContained().size();
1300 return cell_numbers;
std::vector< Box< DIM > > mBoxes
std::vector< unsigned > & rGetHaloNodesRight()
std::set< unsigned > GetLocalBoxes(unsigned boxIndex)
unsigned CalculateContainingBox(Node< DIM > *pNode)
unsigned GetNumLocalBoxes()
~DistributedBoxCollection()
void CalculateInteriorNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs, std::map< unsigned, std::set< unsigned > > &rNodeNeighbours)
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)
void SetCalculateNodeNeighbours(bool calculateNodeNeighbours)
std::set< Node< DIM > * > & rGetNodesContained()
std::map< unsigned, unsigned > mBoxesMapping
std::vector< int > CalculateNumberOfNodesInEachStrip()
unsigned CalculateGlobalIndex(c_vector< unsigned, DIM > coordinateIndices)
unsigned GetProcessOwningNode(Node< DIM > *pNode)
unsigned GetNumLocalRows() const
void SetupLocalBoxesHalfOnly()
c_vector< unsigned, DIM > CalculateCoordinateIndices(unsigned globalIndex)
static const double msFudge
double GetBoxWidth() const
c_vector< unsigned, DIM > mNumBoxesEachDirection
std::vector< unsigned > & rGetHaloNodesLeft()
bool IsOwned(Node< DIM > *pNode)
unsigned GetNumRowsOfBoxes() const
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)
const c_vector< double, SPACE_DIM > & rGetLocation() const
Box< DIM > & rGetHaloBox(unsigned boxIndex)
unsigned mNumBoxesInAFace
Box< DIM > & rGetBox(unsigned boxIndex)
bool GetAreLocalBoxesSet() const
DistributedVectorFactory * mpDistributedBoxStackFactory
void SetupAllLocalBoxes()
void CalculateBoundaryNodePairs(std::vector< Node< DIM > * > &rNodes, std::vector< std::pair< Node< DIM > *, Node< DIM > * > > &rNodePairs, std::map< unsigned, std::set< unsigned > > &rNodeNeighbours)