39 #include "NodePartitioner.hpp"
40 #include "PetscMatTools.hpp"
43 #include "TrianglesMeshReader.hpp"
44 #include "Warnings.hpp"
50 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
52 std::set<unsigned>& rNodesOwned)
59 EXCEPTION(
"The distributed vector factory size in the mesh doesn't match the total number of nodes.");
66 rNodesOwned.insert(node_index);
70 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
72 std::vector<unsigned>& rNodePermutation,
73 std::set<unsigned>& rNodesOwned,
74 std::vector<unsigned>& rProcessorsOffset)
77 assert(ELEMENT_DIM==2 || ELEMENT_DIM==3);
81 WARN_ONCE_ONLY(
"PETSc support for ParMetis is not installed.");
89 unsigned num_local_elements = num_elements / num_procs;
90 unsigned first_local_element = num_local_elements * local_proc_index;
94 num_local_elements += num_elements - (num_local_elements * num_procs);
103 Mat connectivity_matrix;
105 PetscTools::SetupMat(connectivity_matrix, num_nodes, num_nodes, 54, PETSC_DECIDE, PETSC_DECIDE,
false);
110 for (
unsigned element_index = 0; element_index < first_local_element; element_index++)
122 assert(SPACE_DIM >= ELEMENT_DIM);
124 for (
unsigned element_index = 0; element_index < num_local_elements; element_index++)
130 element_data = rMeshReader.
GetElementData(first_local_element + element_index);
137 for (
unsigned i=0; i<element_data.
NodeIndices.size(); i++)
140 for (
unsigned j=0; j<i; j++)
143 MatSetValue(connectivity_matrix, row, col, 1.0, INSERT_VALUES);
144 MatSetValue(connectivity_matrix, col, row, 1.0, INSERT_VALUES);
161 MatGetOwnershipRange(connectivity_matrix, &connectivity_matrix_lo, &connectivity_matrix_hi);
163 unsigned num_local_nodes = connectivity_matrix_hi - connectivity_matrix_lo;
177 MatGetInfo(connectivity_matrix, MAT_LOCAL, &matrix_info);
178 unsigned local_num_nz = (
unsigned) matrix_info.nz_used;
180 size_t size = (num_local_nodes+1)*
sizeof(
PetscInt);
182 PetscMalloc(size, &ptr);
184 size = local_num_nz*
sizeof(
PetscInt);
185 PetscMalloc(size, &ptr);
192 for (
PetscInt row_global_index=connectivity_matrix_lo; row_global_index<connectivity_matrix_hi; row_global_index++)
194 MatGetRow(connectivity_matrix, row_global_index, &row_num_nz, &column_indices, PETSC_NULL);
196 unsigned row_local_index = row_global_index - connectivity_matrix_lo;
197 xadj[row_local_index+1] = xadj[row_local_index] + row_num_nz;
198 for (
PetscInt col_index=0; col_index<row_num_nz; col_index++)
200 adjncy[xadj[row_local_index] + col_index] = column_indices[col_index];
203 MatRestoreRow(connectivity_matrix, row_global_index, &row_num_nz,&column_indices, PETSC_NULL);
210 MatCreateMPIAdj(PETSC_COMM_WORLD, num_local_nodes, num_nodes, xadj, adjncy, PETSC_NULL, &adj_matrix);
219 MatPartitioning part;
220 MatPartitioningCreate(PETSC_COMM_WORLD, &part);
221 MatPartitioningSetAdjacency(part, adj_matrix);
222 MatPartitioningSetFromOptions(part);
223 IS new_process_numbers;
225 MatPartitioningApply(part, &new_process_numbers);
239 #if (PETSC_VERSION_MAJOR == 3) //PETSc 3.x.x
240 ISPartitioningCount(new_process_numbers, num_procs, num_nodes_per_process);
242 ISPartitioningCount(new_process_numbers, num_nodes_per_process);
245 rProcessorsOffset.resize(num_procs);
246 rProcessorsOffset[0] = 0;
247 for (
PetscInt i=1; i<num_procs; i++)
249 rProcessorsOffset[i] = rProcessorsOffset[i-1] + num_nodes_per_process[i-1];
251 unsigned my_num_nodes = num_nodes_per_process[local_proc_index];
252 delete[] num_nodes_per_process;
255 IS new_global_node_indices;
256 ISPartitioningToNumbering(new_process_numbers, &new_global_node_indices);
260 AOCreateBasicIS(new_global_node_indices, PETSC_NULL , &ordering);
264 for (
unsigned i=0; i<my_num_nodes; i++)
266 local_range[i] = rProcessorsOffset[local_proc_index] + i;
268 AOApplicationToPetsc(ordering, my_num_nodes, local_range);
272 rNodesOwned.insert(local_range, local_range + my_num_nodes);
273 delete[] local_range;
277 for (
unsigned i=0; i<num_nodes; i++)
281 AOPetscToApplication(ordering, num_nodes, global_range);
283 rNodePermutation.resize(num_nodes);
284 std::copy(global_range, global_range+num_nodes, rNodePermutation.begin());
285 delete[] global_range;
298 template <
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
300 std::vector<unsigned>& rNodePermutation,
301 std::set<unsigned>& rNodesOwned,
302 std::vector<unsigned>& rProcessorsOffset,
310 std::vector<unsigned> node_ownership(num_nodes, 0);
311 std::vector<unsigned> node_index_ownership(num_nodes,
UNSIGNED_UNSET);
313 for (
unsigned node=0; node < num_nodes; node++)
315 std::vector<double> location = rMeshReader.
GetNextNode();
318 assert(location.size() == SPACE_DIM);
325 bool does_contain =
true;
327 for (
unsigned d=0; d<SPACE_DIM; d++)
330 boundary_check = ((location[d] > lower[d]) || sqrt((location[d]-lower[d])*(location[d]-lower[d])) < DBL_EPSILON);
331 does_contain = (does_contain && boundary_check && (location[d] < upper[d]));
336 node_ownership[node] = 1;
337 rNodesOwned.insert(node);
343 std::vector<unsigned> global_ownership(num_nodes, 0);
344 std::vector<unsigned> global_index_ownership(num_nodes,
UNSIGNED_UNSET);
346 MPI_Allreduce(&node_ownership[0], &global_ownership[0], num_nodes, MPI_UNSIGNED, MPI_LXOR, PETSC_COMM_WORLD);
347 for (
unsigned i=0; i<num_nodes; i++)
349 if (global_ownership[i] == 0)
351 EXCEPTION(
"A node is either not in geometric region, or the regions are not disjoint.");
356 MPI_Allreduce(&node_index_ownership[0], &global_index_ownership[0], num_nodes, MPI_UNSIGNED, MPI_MIN, PETSC_COMM_WORLD);
360 rProcessorsOffset.push_back(rNodePermutation.size());
361 for (
unsigned node=0; node<num_nodes; node++)
363 if (global_index_ownership[node] == proc)
365 rNodePermutation.push_back(node);
369 assert(rNodePermutation.size() == num_nodes);
virtual DistributedVectorFactory * GetDistributedVectorFactory()
virtual ElementData GetNextElementData()=0
virtual ElementData GetElementData(unsigned index)
#define EXCEPTION(message)
static void DumbPartitioning(AbstractMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, std::set< unsigned > &rNodesOwned)
static void GeometricPartitioning(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader, std::vector< unsigned > &rNodePermutation, std::set< unsigned > &rNodesOwned, std::vector< unsigned > &rProcessorsOffset, ChasteCuboid< SPACE_DIM > *pRegion)
virtual unsigned GetNumNodes() const
static void PrintAndReset(std::string message)
std::vector< unsigned > NodeIndices
const ChastePoint< SPACE_DIM > & rGetUpperCorner() const
static void PetscMatrixPartitioning(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader, std::vector< unsigned > &rNodePermutation, std::set< unsigned > &rNodesOwned, std::vector< unsigned > &rProcessorsOffset)
virtual unsigned GetNumElements() const =0
unsigned GetProblemSize() const
const unsigned UNSIGNED_UNSET
virtual bool IsFileFormatBinary()
virtual std::vector< double > GetNextNode()=0
virtual unsigned GetNumNodes() const =0
const ChastePoint< SPACE_DIM > & rGetLowerCorner() const