xref: /petsc/src/dm/partitioner/impls/matpart/partmatpart.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
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 
PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning(PetscPartitioner part,MatPartitioning * mp)8d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning(PetscPartitioner part, MatPartitioning *mp)
9d71ae5a4SJacob Faibussowitsch {
10abe9303eSLisandro Dalcin   PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data;
11abe9303eSLisandro Dalcin 
12abe9303eSLisandro Dalcin   PetscFunctionBegin;
13abe9303eSLisandro Dalcin   *mp = p->mp;
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15abe9303eSLisandro Dalcin }
16abe9303eSLisandro Dalcin 
175d83a8b1SBarry Smith /*@
185d83a8b1SBarry Smith   PetscPartitionerMatPartitioningGetMatPartitioning - Get a `MatPartitioning` instance wrapped by this `PetscPartitioner`.
19abe9303eSLisandro Dalcin 
20abe9303eSLisandro Dalcin   Not Collective
21abe9303eSLisandro Dalcin 
222fe279fdSBarry Smith   Input Parameter:
235d83a8b1SBarry Smith . part - The `PetscPartitioner`
24abe9303eSLisandro Dalcin 
252fe279fdSBarry Smith   Output Parameter:
265d83a8b1SBarry Smith . mp - The `MatPartitioning`
27abe9303eSLisandro Dalcin 
28abe9303eSLisandro Dalcin   Level: developer
29abe9303eSLisandro Dalcin 
30a4e35b19SJacob Faibussowitsch .seealso: `DMPlexDistribute()`, `PetscPartitionerCreate()`
31abe9303eSLisandro Dalcin @*/
PetscPartitionerMatPartitioningGetMatPartitioning(PetscPartitioner part,MatPartitioning * mp)32d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning(PetscPartitioner part, MatPartitioning *mp)
33d71ae5a4SJacob Faibussowitsch {
34abe9303eSLisandro Dalcin   PetscFunctionBegin;
35abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
364f572ea9SToby Isaac   PetscAssertPointer(mp, 2);
37cac4c232SBarry Smith   PetscUseMethod(part, "PetscPartitionerMatPartitioningGetMatPartitioning_C", (PetscPartitioner, MatPartitioning *), (part, mp));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39abe9303eSLisandro Dalcin }
40abe9303eSLisandro Dalcin 
PetscPartitionerDestroy_MatPartitioning(PetscPartitioner part)41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_MatPartitioning(PetscPartitioner part)
42d71ae5a4SJacob Faibussowitsch {
43abe9303eSLisandro Dalcin   PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data;
44abe9303eSLisandro Dalcin 
45abe9303eSLisandro Dalcin   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&p->mp));
472e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)part, "PetscPartitionerMatPartitioningGetMatPartitioning_C", NULL));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50abe9303eSLisandro Dalcin }
51abe9303eSLisandro Dalcin 
PetscPartitionerView_MatPartitioning_ASCII(PetscPartitioner part,PetscViewer viewer)52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_MatPartitioning_ASCII(PetscPartitioner part, PetscViewer viewer)
53d71ae5a4SJacob Faibussowitsch {
54abe9303eSLisandro Dalcin   PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data;
55abe9303eSLisandro Dalcin   PetscViewerFormat                 format;
56abe9303eSLisandro Dalcin 
57abe9303eSLisandro Dalcin   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
599566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MatPartitioning Graph Partitioner:\n"));
609566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
619566063dSJacob Faibussowitsch   if (p->mp) PetscCall(MatPartitioningView(p->mp, viewer));
629566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64abe9303eSLisandro Dalcin }
65abe9303eSLisandro Dalcin 
PetscPartitionerView_MatPartitioning(PetscPartitioner part,PetscViewer viewer)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_MatPartitioning(PetscPartitioner part, PetscViewer viewer)
67d71ae5a4SJacob Faibussowitsch {
68*9f196a02SMartin Diehl   PetscBool isascii;
69abe9303eSLisandro Dalcin 
70abe9303eSLisandro Dalcin   PetscFunctionBegin;
71abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
72abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
73*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
74*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscPartitionerView_MatPartitioning_ASCII(part, viewer));
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76abe9303eSLisandro Dalcin }
77abe9303eSLisandro Dalcin 
PetscPartitionerSetFromOptions_MatPartitioning(PetscPartitioner part,PetscOptionItems PetscOptionsObject)78ce78bad3SBarry Smith static PetscErrorCode PetscPartitionerSetFromOptions_MatPartitioning(PetscPartitioner part, PetscOptionItems PetscOptionsObject)
79d71ae5a4SJacob Faibussowitsch {
80abe9303eSLisandro Dalcin   PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data;
81abe9303eSLisandro Dalcin 
82abe9303eSLisandro Dalcin   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)p->mp, ((PetscObject)part)->prefix));
849566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(p->mp));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86abe9303eSLisandro Dalcin }
87abe9303eSLisandro Dalcin 
PetscPartitionerPartition_MatPartitioning(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertSection,PetscSection edgeSection,PetscSection targetSection,PetscSection partSection,IS * is)8821c92275SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_MatPartitioning(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *is)
89d71ae5a4SJacob Faibussowitsch {
90abe9303eSLisandro Dalcin   PetscPartitioner_MatPartitioning *p = (PetscPartitioner_MatPartitioning *)part->data;
91abe9303eSLisandro Dalcin   Mat                               matadj;
92abe9303eSLisandro Dalcin   IS                                is1, is2, is3;
93abe9303eSLisandro Dalcin   PetscReal                        *tpwgts = NULL;
94abe9303eSLisandro Dalcin   PetscInt                          numVerticesGlobal, numEdges;
95abe9303eSLisandro Dalcin   PetscInt                         *i, *j, *vwgt = NULL;
96abe9303eSLisandro Dalcin   MPI_Comm                          comm;
97abe9303eSLisandro Dalcin 
98abe9303eSLisandro Dalcin   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
100abe9303eSLisandro Dalcin 
101abe9303eSLisandro Dalcin   /* TODO: MatCreateMPIAdj should maybe take global number of ROWS */
102abe9303eSLisandro Dalcin   /* TODO: And vertex distribution in PetscPartitionerPartition_ParMetis should be done using PetscSplitOwnership */
103abe9303eSLisandro Dalcin   numVerticesGlobal = PETSC_DECIDE;
1049566063dSJacob Faibussowitsch   PetscCall(PetscSplitOwnership(comm, &numVertices, &numVerticesGlobal));
105abe9303eSLisandro Dalcin 
106abe9303eSLisandro Dalcin   /* copy arrays to avoid memory errors because MatMPIAdjSetPreallocation copies just pointers */
107abe9303eSLisandro Dalcin   numEdges = start[numVertices];
1089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numVertices + 1, &i));
1099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numEdges, &j));
1109566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(i, start, numVertices + 1));
1119566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(j, adjacency, numEdges));
112abe9303eSLisandro Dalcin 
113abe9303eSLisandro Dalcin   /* construct the adjacency matrix */
1149566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(comm, numVertices, numVerticesGlobal, i, j, NULL, &matadj));
1159566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(p->mp, matadj));
1169566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetNParts(p->mp, nparts));
117abe9303eSLisandro Dalcin 
118abe9303eSLisandro Dalcin   /* calculate partition weights */
119abe9303eSLisandro Dalcin   if (targetSection) {
120abe9303eSLisandro Dalcin     PetscReal sumt;
121abe9303eSLisandro Dalcin     PetscInt  p;
122abe9303eSLisandro Dalcin 
123abe9303eSLisandro Dalcin     sumt = 0.0;
1249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nparts, &tpwgts));
125abe9303eSLisandro Dalcin     for (p = 0; p < nparts; ++p) {
126abe9303eSLisandro Dalcin       PetscInt tpd;
127abe9303eSLisandro Dalcin 
1289566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(targetSection, p, &tpd));
129abe9303eSLisandro Dalcin       sumt += tpd;
130abe9303eSLisandro Dalcin       tpwgts[p] = tpd;
131abe9303eSLisandro Dalcin     }
132abe9303eSLisandro Dalcin     if (sumt) { /* METIS/ParMETIS do not like exactly zero weight */
133abe9303eSLisandro Dalcin       for (p = 0, sumt = 0.0; p < nparts; ++p) {
134abe9303eSLisandro Dalcin         tpwgts[p] = PetscMax(tpwgts[p], PETSC_SMALL);
135abe9303eSLisandro Dalcin         sumt += tpwgts[p];
136abe9303eSLisandro Dalcin       }
137abe9303eSLisandro Dalcin       for (p = 0; p < nparts; ++p) tpwgts[p] /= sumt;
138abe9303eSLisandro Dalcin       for (p = 0, sumt = 0.0; p < nparts - 1; ++p) sumt += tpwgts[p];
139abe9303eSLisandro Dalcin       tpwgts[nparts - 1] = 1. - sumt;
140abe9303eSLisandro Dalcin     } else {
1419566063dSJacob Faibussowitsch       PetscCall(PetscFree(tpwgts));
142abe9303eSLisandro Dalcin     }
143abe9303eSLisandro Dalcin   }
1449566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetPartitionWeights(p->mp, tpwgts));
145abe9303eSLisandro Dalcin 
146abe9303eSLisandro Dalcin   /* calculate vertex weights */
147abe9303eSLisandro Dalcin   if (vertSection) {
148abe9303eSLisandro Dalcin     PetscInt v;
149abe9303eSLisandro Dalcin 
1509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numVertices, &vwgt));
15148a46eb9SPierre Jolivet     for (v = 0; v < numVertices; ++v) PetscCall(PetscSectionGetDof(vertSection, v, &vwgt[v]));
152abe9303eSLisandro Dalcin   }
1539566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetVertexWeights(p->mp, vwgt));
154abe9303eSLisandro Dalcin 
155abe9303eSLisandro Dalcin   /* apply the partitioning */
1569566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(p->mp, &is1));
157abe9303eSLisandro Dalcin 
158abe9303eSLisandro Dalcin   /* construct the PetscSection */
159abe9303eSLisandro Dalcin   {
160abe9303eSLisandro Dalcin     PetscInt        v;
161abe9303eSLisandro Dalcin     const PetscInt *assignment_arr;
162abe9303eSLisandro Dalcin 
1639566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(is1, &assignment_arr));
1649566063dSJacob Faibussowitsch     for (v = 0; v < numVertices; ++v) PetscCall(PetscSectionAddDof(partSection, assignment_arr[v], 1));
1659566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(is1, &assignment_arr));
166abe9303eSLisandro Dalcin   }
167abe9303eSLisandro Dalcin 
168abe9303eSLisandro Dalcin   /* convert assignment IS to global numbering IS */
1699566063dSJacob Faibussowitsch   PetscCall(ISPartitioningToNumbering(is1, &is2));
1709566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
171abe9303eSLisandro Dalcin 
172abe9303eSLisandro Dalcin   /* renumber IS into local numbering */
1739566063dSJacob Faibussowitsch   PetscCall(ISOnComm(is2, PETSC_COMM_SELF, PETSC_USE_POINTER, &is1));
1749566063dSJacob Faibussowitsch   PetscCall(ISRenumber(is1, NULL, NULL, &is3));
1759566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
1769566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
177abe9303eSLisandro Dalcin 
178abe9303eSLisandro Dalcin   /* invert IS */
1799566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(is3));
1809566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(is3, numVertices, &is1));
1819566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is3));
182abe9303eSLisandro Dalcin 
1839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&matadj));
184abe9303eSLisandro Dalcin   *is = is1;
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
186abe9303eSLisandro Dalcin }
187abe9303eSLisandro Dalcin 
PetscPartitionerInitialize_MatPartitioning(PetscPartitioner part)188d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_MatPartitioning(PetscPartitioner part)
189d71ae5a4SJacob Faibussowitsch {
190abe9303eSLisandro Dalcin   PetscFunctionBegin;
191abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_MatPartitioning;
192abe9303eSLisandro Dalcin   part->ops->setfromoptions = PetscPartitionerSetFromOptions_MatPartitioning;
193abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_MatPartitioning;
194abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_MatPartitioning;
1959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)part, "PetscPartitionerMatPartitioningGetMatPartitioning_C", PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning));
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197abe9303eSLisandro Dalcin }
198abe9303eSLisandro Dalcin 
199abe9303eSLisandro Dalcin /*MC
200abe9303eSLisandro Dalcin   PETSCPARTITIONERMATPARTITIONING = "matpartitioning" - A PetscPartitioner object
201abe9303eSLisandro Dalcin 
202abe9303eSLisandro Dalcin   Level: developer
203abe9303eSLisandro Dalcin 
204db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
205abe9303eSLisandro Dalcin M*/
206abe9303eSLisandro Dalcin 
PetscPartitionerCreate_MatPartitioning(PetscPartitioner part)207d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_MatPartitioning(PetscPartitioner part)
208d71ae5a4SJacob Faibussowitsch {
209abe9303eSLisandro Dalcin   PetscPartitioner_MatPartitioning *p;
210abe9303eSLisandro Dalcin 
211abe9303eSLisandro Dalcin   PetscFunctionBegin;
212abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
2134dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
214abe9303eSLisandro Dalcin   part->data = p;
2159566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_MatPartitioning(part));
2169566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(PetscObjectComm((PetscObject)part), &p->mp));
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218abe9303eSLisandro Dalcin }
219