Actual source code: petscpartmatpart.c
 
   petsc-3.10.3 2018-12-18
   
  1:  #include <petsc/private/dmpleximpl.h>
  3: typedef struct {
  4:   MatPartitioning mp;
  5: } PetscPartitioner_MatPartitioning;
  7: static PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning(PetscPartitioner part, MatPartitioning *mp)
  8: {
  9:   PetscPartitioner_MatPartitioning  *p = (PetscPartitioner_MatPartitioning *) part->data;
 12:   *mp = p->mp;
 13:   return(0);
 14: }
 16: /*@C
 17:   PetscPartitionerMatPartitioningGetMatPartitioning - Get a MatPartitioning instance wrapped by this PetscPartitioner.
 19:   Not Collective
 21:   Input Parameters:
 22: . part     - The PetscPartitioner
 24:   Output Parameters:
 25: . mp       - The MatPartitioning
 27:   Level: developer
 29: .seealso PetscPartitionerMatPartitioningSetMatPartitioning(), DMPlexDistribute(), PetscPartitionerCreate()
 30: @*/
 31: PetscErrorCode PetscPartitionerMatPartitioningGetMatPartitioning(PetscPartitioner part, MatPartitioning *mp)
 32: {
 33:   PetscErrorCode          ierr;
 38:   PetscUseMethod(part,"PetscPartitionerMatPartitioningGetMatPartitioning_C",(PetscPartitioner,MatPartitioning*),(part,mp));
 39:   return(0);
 40: }
 42: static PetscErrorCode PetscPartitionerDestroy_MatPartitioning(PetscPartitioner part)
 43: {
 44:   PetscPartitioner_MatPartitioning  *p = (PetscPartitioner_MatPartitioning *) part->data;
 45:   PetscErrorCode                    ierr;
 48:   MatPartitioningDestroy(&p->mp);
 49:   PetscFree(p);
 50:   return(0);
 51: }
 53: static PetscErrorCode PetscPartitionerView_MatPartitioning_Ascii(PetscPartitioner part, PetscViewer viewer)
 54: {
 55:   PetscPartitioner_MatPartitioning  *p = (PetscPartitioner_MatPartitioning *) part->data;
 56:   PetscViewerFormat                 format;
 57:   PetscErrorCode                    ierr;
 60:   PetscViewerGetFormat(viewer, &format);
 61:   PetscViewerASCIIPrintf(viewer, "MatPartitioning Graph Partitioner:\n");
 62:   PetscViewerASCIIPushTab(viewer);
 63:   if (p->mp) {MatPartitioningView(p->mp, viewer);}
 64:   PetscViewerASCIIPopTab(viewer);
 65:   return(0);
 66: }
 68: static PetscErrorCode PetscPartitionerView_MatPartitioning(PetscPartitioner part, PetscViewer viewer)
 69: {
 70:   PetscBool      iascii;
 76:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 77:   if (iascii) {PetscPartitionerView_MatPartitioning_Ascii(part, viewer);}
 78:   return(0);
 79: }
 81: static PetscErrorCode PetscPartitionerSetFromOptions_MatPartitioning(PetscOptionItems *PetscOptionsObject, PetscPartitioner part)
 82: {
 83:   PetscPartitioner_MatPartitioning  *p = (PetscPartitioner_MatPartitioning *) part->data;
 84:   PetscErrorCode                    ierr;
 87:   PetscObjectSetOptionsPrefix((PetscObject)p->mp,((PetscObject)part)->prefix);
 88:   MatPartitioningSetFromOptions(p->mp);
 89:   return(0);
 90: }
 92: static PetscErrorCode PetscPartitionerPartition_MatPartitioning(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *is)
 93: {
 94:   PetscPartitioner_MatPartitioning  *p = (PetscPartitioner_MatPartitioning *) part->data;
 95:   Mat                               matadj;
 96:   IS                                is1, is2, is3;
 97:   PetscInt                          numVerticesGlobal, numEdges;
 98:   PetscInt                          *i, *j;
 99:   MPI_Comm                          comm;
100:   PetscErrorCode                    ierr;
103:   PetscObjectGetComm((PetscObject)part, &comm);
104:   if (numVertices < 0) SETERRQ(comm, PETSC_ERR_PLIB, "number of vertices must be specified");
106:   /* TODO: MatCreateMPIAdj should maybe take global number of ROWS */
107:   /* TODO: And vertex distribution in PetscPartitionerPartition_ParMetis should be done using PetscSplitOwnership */
108:   numVerticesGlobal = PETSC_DECIDE;
109:   PetscSplitOwnership(comm, &numVertices, &numVerticesGlobal);
111:   /* copy arrays to avoid memory errors because MatMPIAdjSetPreallocation copies just pointers */
112:   numEdges = start[numVertices];
113:   PetscMalloc1(numVertices+1, &i);
114:   PetscMalloc1(numEdges, &j);
115:   PetscMemcpy(i, start, (numVertices+1)*sizeof(PetscInt));
116:   PetscMemcpy(j, adjacency, numEdges*sizeof(PetscInt));
118:   /* construct the adjacency matrix */
119:   MatCreateMPIAdj(comm, numVertices, numVerticesGlobal, i, j, NULL, &matadj);
120:   MatPartitioningSetAdjacency(p->mp, matadj);
121:   MatPartitioningSetNParts(p->mp, nparts);
123:   /* apply the partitioning */
124:   MatPartitioningApply(p->mp, &is1);
126:   /* construct the PetscSection */
127:   {
128:     PetscInt v;
129:     const PetscInt *assignment_arr;
131:     PetscSectionSetChart(partSection, 0, nparts);
132:     ISGetIndices(is1, &assignment_arr);
133:     for (v = 0; v < numVertices; ++v) {PetscSectionAddDof(partSection, assignment_arr[v], 1);}
134:     ISRestoreIndices(is1, &assignment_arr);
135:     PetscSectionSetUp(partSection);
136:   }
138:   /* convert assignment IS to global numbering IS */
139:   ISPartitioningToNumbering(is1, &is2);
140:   ISDestroy(&is1);
142:   /* renumber IS into local numbering */
143:   ISOnComm(is2, PETSC_COMM_SELF, PETSC_USE_POINTER, &is1);
144:   ISRenumber(is1, NULL, NULL, &is3);
145:   ISDestroy(&is1);
146:   ISDestroy(&is2);
148:   /* invert IS */
149:   ISSetPermutation(is3);
150:   ISInvertPermutation(is3, numVertices, &is1);
151:   ISDestroy(&is3);
153:   MatDestroy(&matadj);
154:   *is = is1;
155:   return(0);
156: }
158: static PetscErrorCode PetscPartitionerInitialize_MatPartitioning(PetscPartitioner part)
159: {
160:   PetscErrorCode                    ierr;
163:   part->ops->view           = PetscPartitionerView_MatPartitioning;
164:   part->ops->setfromoptions = PetscPartitionerSetFromOptions_MatPartitioning;
165:   part->ops->destroy        = PetscPartitionerDestroy_MatPartitioning;
166:   part->ops->partition      = PetscPartitionerPartition_MatPartitioning;
167:   PetscObjectComposeFunction((PetscObject)part,"PetscPartitionerMatPartitioningGetMatPartitioning_C",PetscPartitionerMatPartitioningGetMatPartitioning_MatPartitioning);
168:   return(0);
169: }
171: /*MC
172:   PETSCPARTITIONERMATPARTITIONING = "matpartitioning" - A PetscPartitioner object
174:   Level: developer
176: .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
177: M*/
179: PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_MatPartitioning(PetscPartitioner part)
180: {
181:   PetscPartitioner_MatPartitioning  *p;
182:   PetscErrorCode                    ierr;
186:   PetscNewLog(part, &p);
187:   part->data = p;
188:   PetscPartitionerInitialize_MatPartitioning(part);
189:   MatPartitioningCreate(PetscObjectComm((PetscObject)part), &p->mp);
190:   return(0);
191: }