1abe9303eSLisandro Dalcin #include <petscmat.h> 2abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/ 3abe9303eSLisandro Dalcin 4abe9303eSLisandro Dalcin typedef struct { 5abe9303eSLisandro Dalcin MatPartitioning mp; 6abe9303eSLisandro Dalcin } PetscPartitioner_MatPartitioning; 7abe9303eSLisandro Dalcin 89371c9d4SSatish Balay static PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning(PetscPartitioner part, MatPartitioning *mp) { 9abe9303eSLisandro Dalcin PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data; 10abe9303eSLisandro Dalcin 11abe9303eSLisandro Dalcin PetscFunctionBegin; 12abe9303eSLisandro Dalcin *mp = p->mp; 13abe9303eSLisandro Dalcin PetscFunctionReturn(0); 14abe9303eSLisandro Dalcin } 15abe9303eSLisandro Dalcin 16abe9303eSLisandro Dalcin /*@C 17abe9303eSLisandro Dalcin PetscPartitionerMatPartitioningGetMatPartitioning - Get a MatPartitioning instance wrapped by this PetscPartitioner. 18abe9303eSLisandro Dalcin 19abe9303eSLisandro Dalcin Not Collective 20abe9303eSLisandro Dalcin 21abe9303eSLisandro Dalcin Input Parameters: 22abe9303eSLisandro Dalcin . part - The PetscPartitioner 23abe9303eSLisandro Dalcin 24abe9303eSLisandro Dalcin Output Parameters: 25abe9303eSLisandro Dalcin . mp - The MatPartitioning 26abe9303eSLisandro Dalcin 27abe9303eSLisandro Dalcin Level: developer 28abe9303eSLisandro Dalcin 29db781477SPatrick Sanan .seealso `DMPlexDistribute()`, `PetscPartitionerCreate()` 30abe9303eSLisandro Dalcin @*/ 319371c9d4SSatish Balay PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning(PetscPartitioner part, MatPartitioning *mp) { 32abe9303eSLisandro Dalcin PetscFunctionBegin; 33abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 34abe9303eSLisandro Dalcin PetscValidPointer(mp, 2); 35cac4c232SBarry Smith PetscUseMethod(part, "PetscPartitionerMatPartitioningGetMatPartitioning_C", (PetscPartitioner, MatPartitioning *), (part, mp)); 36abe9303eSLisandro Dalcin PetscFunctionReturn(0); 37abe9303eSLisandro Dalcin } 38abe9303eSLisandro Dalcin 399371c9d4SSatish Balay static PetscErrorCode PetscPartitionerDestroy_MatPartitioning(PetscPartitioner part) { 40abe9303eSLisandro Dalcin PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data; 41abe9303eSLisandro Dalcin 42abe9303eSLisandro Dalcin PetscFunctionBegin; 439566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&p->mp)); 442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)part, "PetscPartitionerMatPartitioningGetMatPartitioning_C", NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data)); 46abe9303eSLisandro Dalcin PetscFunctionReturn(0); 47abe9303eSLisandro Dalcin } 48abe9303eSLisandro Dalcin 499371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_MatPartitioning_ASCII(PetscPartitioner part, PetscViewer viewer) { 50abe9303eSLisandro Dalcin PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data; 51abe9303eSLisandro Dalcin PetscViewerFormat format; 52abe9303eSLisandro Dalcin 53abe9303eSLisandro Dalcin PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MatPartitioning Graph Partitioner:\n")); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 579566063dSJacob Faibussowitsch if (p->mp) PetscCall(MatPartitioningView(p->mp, viewer)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 59abe9303eSLisandro Dalcin PetscFunctionReturn(0); 60abe9303eSLisandro Dalcin } 61abe9303eSLisandro Dalcin 629371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_MatPartitioning(PetscPartitioner part, PetscViewer viewer) { 63abe9303eSLisandro Dalcin PetscBool iascii; 64abe9303eSLisandro Dalcin 65abe9303eSLisandro Dalcin PetscFunctionBegin; 66abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 67abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 699566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscPartitionerView_MatPartitioning_ASCII(part, viewer)); 70abe9303eSLisandro Dalcin PetscFunctionReturn(0); 71abe9303eSLisandro Dalcin } 72abe9303eSLisandro Dalcin 739371c9d4SSatish Balay static PetscErrorCode PetscPartitionerSetFromOptions_MatPartitioning(PetscPartitioner part, PetscOptionItems *PetscOptionsObject) { 74abe9303eSLisandro Dalcin PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data; 75abe9303eSLisandro Dalcin 76abe9303eSLisandro Dalcin PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)p->mp, ((PetscObject)part)->prefix)); 789566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(p->mp)); 79abe9303eSLisandro Dalcin PetscFunctionReturn(0); 80abe9303eSLisandro Dalcin } 81abe9303eSLisandro Dalcin 829371c9d4SSatish Balay static PetscErrorCode PetscPartitionerPartition_MatPartitioning(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *is) { 83abe9303eSLisandro Dalcin PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data; 84abe9303eSLisandro Dalcin Mat matadj; 85abe9303eSLisandro Dalcin IS is1, is2, is3; 86abe9303eSLisandro Dalcin PetscReal *tpwgts = NULL; 87abe9303eSLisandro Dalcin PetscInt numVerticesGlobal, numEdges; 88abe9303eSLisandro Dalcin PetscInt *i, *j, *vwgt = NULL; 89abe9303eSLisandro Dalcin MPI_Comm comm; 90abe9303eSLisandro Dalcin 91abe9303eSLisandro Dalcin PetscFunctionBegin; 929566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)part, &comm)); 93abe9303eSLisandro Dalcin 94abe9303eSLisandro Dalcin /* TODO: MatCreateMPIAdj should maybe take global number of ROWS */ 95abe9303eSLisandro Dalcin /* TODO: And vertex distribution in PetscPartitionerPartition_ParMetis should be done using PetscSplitOwnership */ 96abe9303eSLisandro Dalcin numVerticesGlobal = PETSC_DECIDE; 979566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnership(comm, &numVertices, &numVerticesGlobal)); 98abe9303eSLisandro Dalcin 99abe9303eSLisandro Dalcin /* copy arrays to avoid memory errors because MatMPIAdjSetPreallocation copies just pointers */ 100abe9303eSLisandro Dalcin numEdges = start[numVertices]; 1019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numVertices + 1, &i)); 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numEdges, &j)); 1039566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(i, start, numVertices + 1)); 1049566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(j, adjacency, numEdges)); 105abe9303eSLisandro Dalcin 106abe9303eSLisandro Dalcin /* construct the adjacency matrix */ 1079566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(comm, numVertices, numVerticesGlobal, i, j, NULL, &matadj)); 1089566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(p->mp, matadj)); 1099566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(p->mp, nparts)); 110abe9303eSLisandro Dalcin 111abe9303eSLisandro Dalcin /* calculate partition weights */ 112abe9303eSLisandro Dalcin if (targetSection) { 113abe9303eSLisandro Dalcin PetscReal sumt; 114abe9303eSLisandro Dalcin PetscInt p; 115abe9303eSLisandro Dalcin 116abe9303eSLisandro Dalcin sumt = 0.0; 1179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nparts, &tpwgts)); 118abe9303eSLisandro Dalcin for (p = 0; p < nparts; ++p) { 119abe9303eSLisandro Dalcin PetscInt tpd; 120abe9303eSLisandro Dalcin 1219566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(targetSection, p, &tpd)); 122abe9303eSLisandro Dalcin sumt += tpd; 123abe9303eSLisandro Dalcin tpwgts[p] = tpd; 124abe9303eSLisandro Dalcin } 125abe9303eSLisandro Dalcin if (sumt) { /* METIS/ParMETIS do not like exactly zero weight */ 126abe9303eSLisandro Dalcin for (p = 0, sumt = 0.0; p < nparts; ++p) { 127abe9303eSLisandro Dalcin tpwgts[p] = PetscMax(tpwgts[p], PETSC_SMALL); 128abe9303eSLisandro Dalcin sumt += tpwgts[p]; 129abe9303eSLisandro Dalcin } 130abe9303eSLisandro Dalcin for (p = 0; p < nparts; ++p) tpwgts[p] /= sumt; 131abe9303eSLisandro Dalcin for (p = 0, sumt = 0.0; p < nparts - 1; ++p) sumt += tpwgts[p]; 132abe9303eSLisandro Dalcin tpwgts[nparts - 1] = 1. - sumt; 133abe9303eSLisandro Dalcin } else { 1349566063dSJacob Faibussowitsch PetscCall(PetscFree(tpwgts)); 135abe9303eSLisandro Dalcin } 136abe9303eSLisandro Dalcin } 1379566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetPartitionWeights(p->mp, tpwgts)); 138abe9303eSLisandro Dalcin 139abe9303eSLisandro Dalcin /* calculate vertex weights */ 140abe9303eSLisandro Dalcin if (vertSection) { 141abe9303eSLisandro Dalcin PetscInt v; 142abe9303eSLisandro Dalcin 1439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numVertices, &vwgt)); 144*48a46eb9SPierre Jolivet for (v = 0; v < numVertices; ++v) PetscCall(PetscSectionGetDof(vertSection, v, &vwgt[v])); 145abe9303eSLisandro Dalcin } 1469566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetVertexWeights(p->mp, vwgt)); 147abe9303eSLisandro Dalcin 148abe9303eSLisandro Dalcin /* apply the partitioning */ 1499566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(p->mp, &is1)); 150abe9303eSLisandro Dalcin 151abe9303eSLisandro Dalcin /* construct the PetscSection */ 152abe9303eSLisandro Dalcin { 153abe9303eSLisandro Dalcin PetscInt v; 154abe9303eSLisandro Dalcin const PetscInt *assignment_arr; 155abe9303eSLisandro Dalcin 1569566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is1, &assignment_arr)); 1579566063dSJacob Faibussowitsch for (v = 0; v < numVertices; ++v) PetscCall(PetscSectionAddDof(partSection, assignment_arr[v], 1)); 1589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is1, &assignment_arr)); 159abe9303eSLisandro Dalcin } 160abe9303eSLisandro Dalcin 161abe9303eSLisandro Dalcin /* convert assignment IS to global numbering IS */ 1629566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(is1, &is2)); 1639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 164abe9303eSLisandro Dalcin 165abe9303eSLisandro Dalcin /* renumber IS into local numbering */ 1669566063dSJacob Faibussowitsch PetscCall(ISOnComm(is2, PETSC_COMM_SELF, PETSC_USE_POINTER, &is1)); 1679566063dSJacob Faibussowitsch PetscCall(ISRenumber(is1, NULL, NULL, &is3)); 1689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 1699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 170abe9303eSLisandro Dalcin 171abe9303eSLisandro Dalcin /* invert IS */ 1729566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(is3)); 1739566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(is3, numVertices, &is1)); 1749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is3)); 175abe9303eSLisandro Dalcin 1769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&matadj)); 177abe9303eSLisandro Dalcin *is = is1; 178abe9303eSLisandro Dalcin PetscFunctionReturn(0); 179abe9303eSLisandro Dalcin } 180abe9303eSLisandro Dalcin 1819371c9d4SSatish Balay static PetscErrorCode PetscPartitionerInitialize_MatPartitioning(PetscPartitioner part) { 182abe9303eSLisandro Dalcin PetscFunctionBegin; 183abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_MatPartitioning; 184abe9303eSLisandro Dalcin part->ops->setfromoptions = PetscPartitionerSetFromOptions_MatPartitioning; 185abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_MatPartitioning; 186abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_MatPartitioning; 1879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)part, "PetscPartitionerMatPartitioningGetMatPartitioning_C", PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning)); 188abe9303eSLisandro Dalcin PetscFunctionReturn(0); 189abe9303eSLisandro Dalcin } 190abe9303eSLisandro Dalcin 191abe9303eSLisandro Dalcin /*MC 192abe9303eSLisandro Dalcin PETSCPARTITIONERMATPARTITIONING = "matpartitioning" - A PetscPartitioner object 193abe9303eSLisandro Dalcin 194abe9303eSLisandro Dalcin Level: developer 195abe9303eSLisandro Dalcin 196db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()` 197abe9303eSLisandro Dalcin M*/ 198abe9303eSLisandro Dalcin 1999371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_MatPartitioning(PetscPartitioner part) { 200abe9303eSLisandro Dalcin PetscPartitioner_MatPartitioning *p; 201abe9303eSLisandro Dalcin 202abe9303eSLisandro Dalcin PetscFunctionBegin; 203abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 2049566063dSJacob Faibussowitsch PetscCall(PetscNewLog(part, &p)); 205abe9303eSLisandro Dalcin part->data = p; 2069566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_MatPartitioning(part)); 2079566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PetscObjectComm((PetscObject)part), &p->mp)); 208abe9303eSLisandro Dalcin PetscFunctionReturn(0); 209abe9303eSLisandro Dalcin } 210