1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 270034214SMatthew G. Knepley 377623264SMatthew G. Knepley PetscClassId PETSCPARTITIONER_CLASSID = 0; 477623264SMatthew G. Knepley 577623264SMatthew G. Knepley PetscFunctionList PetscPartitionerList = NULL; 677623264SMatthew G. Knepley PetscBool PetscPartitionerRegisterAllCalled = PETSC_FALSE; 777623264SMatthew G. Knepley 877623264SMatthew G. Knepley PetscBool ChacoPartitionercite = PETSC_FALSE; 977623264SMatthew G. Knepley const char ChacoPartitionerCitation[] = "@inproceedings{Chaco95,\n" 1077623264SMatthew G. Knepley " author = {Bruce Hendrickson and Robert Leland},\n" 1177623264SMatthew G. Knepley " title = {A multilevel algorithm for partitioning graphs},\n" 1277623264SMatthew G. Knepley " booktitle = {Supercomputing '95: Proceedings of the 1995 ACM/IEEE Conference on Supercomputing (CDROM)}," 1377623264SMatthew G. Knepley " isbn = {0-89791-816-9},\n" 1477623264SMatthew G. Knepley " pages = {28},\n" 1577623264SMatthew G. Knepley " doi = {http://doi.acm.org/10.1145/224170.224228},\n" 1677623264SMatthew G. Knepley " publisher = {ACM Press},\n" 1777623264SMatthew G. Knepley " address = {New York},\n" 1877623264SMatthew G. Knepley " year = {1995}\n}\n"; 1977623264SMatthew G. Knepley 2077623264SMatthew G. Knepley PetscBool ParMetisPartitionercite = PETSC_FALSE; 2177623264SMatthew G. Knepley const char ParMetisPartitionerCitation[] = "@article{KarypisKumar98,\n" 2277623264SMatthew G. Knepley " author = {George Karypis and Vipin Kumar},\n" 2377623264SMatthew G. Knepley " title = {A Parallel Algorithm for Multilevel Graph Partitioning and Sparse Matrix Ordering},\n" 2477623264SMatthew G. Knepley " journal = {Journal of Parallel and Distributed Computing},\n" 2577623264SMatthew G. Knepley " volume = {48},\n" 2677623264SMatthew G. Knepley " pages = {71--85},\n" 2777623264SMatthew G. Knepley " year = {1998}\n}\n"; 2877623264SMatthew G. Knepley 2970034214SMatthew G. Knepley #undef __FUNCT__ 30532c4e7dSMichael Lange #define __FUNCT__ "DMPlexCreatePartitionerGraph" 31532c4e7dSMichael Lange /*@C 32532c4e7dSMichael Lange DMPlexCreatePartitionerGraph - Create a CSR graph of point connections for the partitioner 33532c4e7dSMichael Lange 34532c4e7dSMichael Lange Input Parameters: 35532c4e7dSMichael Lange + dm - The mesh DM dm 36532c4e7dSMichael Lange - height - Height of the strata from which to construct the graph 37532c4e7dSMichael Lange 38532c4e7dSMichael Lange Output Parameter: 39532c4e7dSMichael Lange + numVertices - Number of vertices in the graph 40532c4e7dSMichael Lange - offsets - Point offsets in the graph 41532c4e7dSMichael Lange - adjacency - Point connectivity in the graph 42532c4e7dSMichael Lange 43532c4e7dSMichael Lange The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 44532c4e7dSMichael Lange DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 45532c4e7dSMichael Lange representation on the mesh. 46532c4e7dSMichael Lange 47532c4e7dSMichael Lange Level: developer 48532c4e7dSMichael Lange 49532c4e7dSMichael Lange .seealso: PetscPartitionerGetType(), PetscPartitionerCreate(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 50532c4e7dSMichael Lange @*/ 51532c4e7dSMichael Lange PetscErrorCode DMPlexCreatePartitionerGraph(DM dm, PetscInt height, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency) 52532c4e7dSMichael Lange { 5343f7d02bSMichael Lange PetscInt p, pStart, pEnd, a, adjSize, idx, size, nroots; 54389e55d8SMichael Lange PetscInt *adj = NULL, *vOffsets = NULL, *graph = NULL; 558cfe4c1fSMichael Lange IS cellNumbering; 568cfe4c1fSMichael Lange const PetscInt *cellNum; 57532c4e7dSMichael Lange PetscBool useCone, useClosure; 58532c4e7dSMichael Lange PetscSection section; 59532c4e7dSMichael Lange PetscSegBuffer adjBuffer; 608cfe4c1fSMichael Lange PetscSF sfPoint; 61532c4e7dSMichael Lange PetscMPIInt rank; 62532c4e7dSMichael Lange PetscErrorCode ierr; 63532c4e7dSMichael Lange 64532c4e7dSMichael Lange PetscFunctionBegin; 65532c4e7dSMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 66532c4e7dSMichael Lange ierr = DMPlexGetHeightStratum(dm, height, &pStart, &pEnd);CHKERRQ(ierr); 678cfe4c1fSMichael Lange ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 688cfe4c1fSMichael Lange ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 69532c4e7dSMichael Lange /* Build adjacency graph via a section/segbuffer */ 70532c4e7dSMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 71532c4e7dSMichael Lange ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 72532c4e7dSMichael Lange ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&adjBuffer);CHKERRQ(ierr); 73532c4e7dSMichael Lange /* Always use FVM adjacency to create partitioner graph */ 74532c4e7dSMichael Lange ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr); 75532c4e7dSMichael Lange ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr); 76532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr); 77532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);CHKERRQ(ierr); 788cfe4c1fSMichael Lange if (nroots > 0) { 798cfe4c1fSMichael Lange ierr = DMPlexGetCellNumbering(dm, &cellNumbering);CHKERRQ(ierr); 808cfe4c1fSMichael Lange ierr = ISGetIndices(cellNumbering, &cellNum);CHKERRQ(ierr); 818cfe4c1fSMichael Lange } 828cfe4c1fSMichael Lange for (*numVertices = 0, p = pStart; p < pEnd; p++) { 838cfe4c1fSMichael Lange /* Skip non-owned cells in parallel (ParMetis expects no overlap) */ 848cfe4c1fSMichael Lange if (nroots > 0) {if (cellNum[p] < 0) continue;} 85532c4e7dSMichael Lange adjSize = PETSC_DETERMINE; 86532c4e7dSMichael Lange ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 87532c4e7dSMichael Lange for (a = 0; a < adjSize; ++a) { 88532c4e7dSMichael Lange const PetscInt point = adj[a]; 89532c4e7dSMichael Lange if (point != p && pStart <= point && point < pEnd) { 90532c4e7dSMichael Lange PetscInt *PETSC_RESTRICT pBuf; 91532c4e7dSMichael Lange ierr = PetscSectionAddDof(section, p, 1);CHKERRQ(ierr); 92532c4e7dSMichael Lange ierr = PetscSegBufferGetInts(adjBuffer, 1, &pBuf);CHKERRQ(ierr); 93532c4e7dSMichael Lange *pBuf = point; 94532c4e7dSMichael Lange } 95532c4e7dSMichael Lange } 968cfe4c1fSMichael Lange (*numVertices)++; 97532c4e7dSMichael Lange } 98532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr); 99532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseClosure(dm, useClosure);CHKERRQ(ierr); 100532c4e7dSMichael Lange /* Derive CSR graph from section/segbuffer */ 101532c4e7dSMichael Lange ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 102532c4e7dSMichael Lange ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 103389e55d8SMichael Lange ierr = PetscMalloc1(*numVertices+1, &vOffsets);CHKERRQ(ierr); 10443f7d02bSMichael Lange for (idx = 0, p = pStart; p < pEnd; p++) { 10543f7d02bSMichael Lange if (nroots > 0) {if (cellNum[p] < 0) continue;} 10643f7d02bSMichael Lange ierr = PetscSectionGetOffset(section, p, &(vOffsets[idx++]));CHKERRQ(ierr); 10743f7d02bSMichael Lange } 108532c4e7dSMichael Lange vOffsets[*numVertices] = size; 109532c4e7dSMichael Lange if (offsets) *offsets = vOffsets; 110389e55d8SMichael Lange ierr = PetscSegBufferExtractAlloc(adjBuffer, &graph);CHKERRQ(ierr); 1118cfe4c1fSMichael Lange if (nroots > 0) { 1128cfe4c1fSMichael Lange ISLocalToGlobalMapping ltogCells; 1138cfe4c1fSMichael Lange PetscInt n, size, *cells_arr; 1148cfe4c1fSMichael Lange /* In parallel, apply a global cell numbering to the graph */ 1158cfe4c1fSMichael Lange ierr = ISRestoreIndices(cellNumbering, &cellNum);CHKERRQ(ierr); 1168cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingCreateIS(cellNumbering, <ogCells);CHKERRQ(ierr); 1178cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingGetSize(ltogCells, &size);CHKERRQ(ierr); 1188cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingGetIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr); 1198cfe4c1fSMichael Lange /* Convert to positive global cell numbers */ 1208cfe4c1fSMichael Lange for (n=0; n<size; n++) {if (cells_arr[n] < 0) cells_arr[n] = -(cells_arr[n]+1);} 1218cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingRestoreIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr); 122389e55d8SMichael Lange ierr = ISLocalToGlobalMappingApplyBlock(ltogCells, vOffsets[*numVertices], graph, graph);CHKERRQ(ierr); 1238cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingDestroy(<ogCells);CHKERRQ(ierr); 1248cfe4c1fSMichael Lange } 125389e55d8SMichael Lange if (adjacency) *adjacency = graph; 126532c4e7dSMichael Lange /* Clean up */ 127532c4e7dSMichael Lange ierr = PetscSegBufferDestroy(&adjBuffer);CHKERRQ(ierr); 128532c4e7dSMichael Lange ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 129532c4e7dSMichael Lange ierr = PetscFree(adj);CHKERRQ(ierr); 130532c4e7dSMichael Lange PetscFunctionReturn(0); 131532c4e7dSMichael Lange } 132532c4e7dSMichael Lange 133532c4e7dSMichael Lange #undef __FUNCT__ 13470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateNeighborCSR" 13570034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency) 13670034214SMatthew G. Knepley { 13770034214SMatthew G. Knepley const PetscInt maxFaceCases = 30; 13870034214SMatthew G. Knepley PetscInt numFaceCases = 0; 13970034214SMatthew G. Knepley PetscInt numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */ 14070034214SMatthew G. Knepley PetscInt *off, *adj; 14170034214SMatthew G. Knepley PetscInt *neighborCells = NULL; 14270034214SMatthew G. Knepley PetscInt dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell; 14370034214SMatthew G. Knepley PetscErrorCode ierr; 14470034214SMatthew G. Knepley 14570034214SMatthew G. Knepley PetscFunctionBegin; 14670034214SMatthew G. Knepley /* For parallel partitioning, I think you have to communicate supports */ 147c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 14870034214SMatthew G. Knepley cellDim = dim - cellHeight; 14970034214SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 15070034214SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 15170034214SMatthew G. Knepley if (cEnd - cStart == 0) { 15270034214SMatthew G. Knepley if (numVertices) *numVertices = 0; 15370034214SMatthew G. Knepley if (offsets) *offsets = NULL; 15470034214SMatthew G. Knepley if (adjacency) *adjacency = NULL; 15570034214SMatthew G. Knepley PetscFunctionReturn(0); 15670034214SMatthew G. Knepley } 15770034214SMatthew G. Knepley numCells = cEnd - cStart; 15870034214SMatthew G. Knepley faceDepth = depth - cellHeight; 15970034214SMatthew G. Knepley if (dim == depth) { 16070034214SMatthew G. Knepley PetscInt f, fStart, fEnd; 16170034214SMatthew G. Knepley 16270034214SMatthew G. Knepley ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr); 16370034214SMatthew G. Knepley /* Count neighboring cells */ 16470034214SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr); 16570034214SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 16670034214SMatthew G. Knepley const PetscInt *support; 16770034214SMatthew G. Knepley PetscInt supportSize; 16870034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 16970034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr); 17070034214SMatthew G. Knepley if (supportSize == 2) { 1719ffc88e4SToby Isaac PetscInt numChildren; 1729ffc88e4SToby Isaac 1739ffc88e4SToby Isaac ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 1749ffc88e4SToby Isaac if (!numChildren) { 17570034214SMatthew G. Knepley ++off[support[0]-cStart+1]; 17670034214SMatthew G. Knepley ++off[support[1]-cStart+1]; 17770034214SMatthew G. Knepley } 17870034214SMatthew G. Knepley } 1799ffc88e4SToby Isaac } 18070034214SMatthew G. Knepley /* Prefix sum */ 18170034214SMatthew G. Knepley for (c = 1; c <= numCells; ++c) off[c] += off[c-1]; 18270034214SMatthew G. Knepley if (adjacency) { 18370034214SMatthew G. Knepley PetscInt *tmp; 18470034214SMatthew G. Knepley 18570034214SMatthew G. Knepley ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr); 186854ce69bSBarry Smith ierr = PetscMalloc1(numCells+1, &tmp);CHKERRQ(ierr); 18770034214SMatthew G. Knepley ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr); 18870034214SMatthew G. Knepley /* Get neighboring cells */ 18970034214SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 19070034214SMatthew G. Knepley const PetscInt *support; 19170034214SMatthew G. Knepley PetscInt supportSize; 19270034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 19370034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr); 19470034214SMatthew G. Knepley if (supportSize == 2) { 1959ffc88e4SToby Isaac PetscInt numChildren; 1969ffc88e4SToby Isaac 1979ffc88e4SToby Isaac ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 1989ffc88e4SToby Isaac if (!numChildren) { 19970034214SMatthew G. Knepley adj[tmp[support[0]-cStart]++] = support[1]; 20070034214SMatthew G. Knepley adj[tmp[support[1]-cStart]++] = support[0]; 20170034214SMatthew G. Knepley } 20270034214SMatthew G. Knepley } 2039ffc88e4SToby Isaac } 20470034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG) 20570034214SMatthew G. Knepley for (c = 0; c < cEnd-cStart; ++c) if (tmp[c] != off[c+1]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Offset %d != %d for cell %d", tmp[c], off[c], c+cStart); 20670034214SMatthew G. Knepley #endif 20770034214SMatthew G. Knepley ierr = PetscFree(tmp);CHKERRQ(ierr); 20870034214SMatthew G. Knepley } 20970034214SMatthew G. Knepley if (numVertices) *numVertices = numCells; 21070034214SMatthew G. Knepley if (offsets) *offsets = off; 21170034214SMatthew G. Knepley if (adjacency) *adjacency = adj; 21270034214SMatthew G. Knepley PetscFunctionReturn(0); 21370034214SMatthew G. Knepley } 21470034214SMatthew G. Knepley /* Setup face recognition */ 21570034214SMatthew G. Knepley if (faceDepth == 1) { 21670034214SMatthew G. Knepley PetscInt cornersSeen[30] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0}; /* Could use PetscBT */ 21770034214SMatthew G. Knepley 21870034214SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 21970034214SMatthew G. Knepley PetscInt corners; 22070034214SMatthew G. Knepley 22170034214SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr); 22270034214SMatthew G. Knepley if (!cornersSeen[corners]) { 22370034214SMatthew G. Knepley PetscInt nFV; 22470034214SMatthew G. Knepley 22570034214SMatthew G. Knepley if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases"); 22670034214SMatthew G. Knepley cornersSeen[corners] = 1; 22770034214SMatthew G. Knepley 22870034214SMatthew G. Knepley ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr); 22970034214SMatthew G. Knepley 23070034214SMatthew G. Knepley numFaceVertices[numFaceCases++] = nFV; 23170034214SMatthew G. Knepley } 23270034214SMatthew G. Knepley } 23370034214SMatthew G. Knepley } 23470034214SMatthew G. Knepley ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr); 23570034214SMatthew G. Knepley /* Count neighboring cells */ 23670034214SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 23770034214SMatthew G. Knepley PetscInt numNeighbors = PETSC_DETERMINE, n; 23870034214SMatthew G. Knepley 2398b0b4c70SToby Isaac ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr); 24070034214SMatthew G. Knepley /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */ 24170034214SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 24270034214SMatthew G. Knepley PetscInt cellPair[2]; 24370034214SMatthew G. Knepley PetscBool found = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE; 24470034214SMatthew G. Knepley PetscInt meetSize = 0; 24570034214SMatthew G. Knepley const PetscInt *meet = NULL; 24670034214SMatthew G. Knepley 24770034214SMatthew G. Knepley cellPair[0] = cell; cellPair[1] = neighborCells[n]; 24870034214SMatthew G. Knepley if (cellPair[0] == cellPair[1]) continue; 24970034214SMatthew G. Knepley if (!found) { 25070034214SMatthew G. Knepley ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 25170034214SMatthew G. Knepley if (meetSize) { 25270034214SMatthew G. Knepley PetscInt f; 25370034214SMatthew G. Knepley 25470034214SMatthew G. Knepley for (f = 0; f < numFaceCases; ++f) { 25570034214SMatthew G. Knepley if (numFaceVertices[f] == meetSize) { 25670034214SMatthew G. Knepley found = PETSC_TRUE; 25770034214SMatthew G. Knepley break; 25870034214SMatthew G. Knepley } 25970034214SMatthew G. Knepley } 26070034214SMatthew G. Knepley } 26170034214SMatthew G. Knepley ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 26270034214SMatthew G. Knepley } 26370034214SMatthew G. Knepley if (found) ++off[cell-cStart+1]; 26470034214SMatthew G. Knepley } 26570034214SMatthew G. Knepley } 26670034214SMatthew G. Knepley /* Prefix sum */ 26770034214SMatthew G. Knepley for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1]; 26870034214SMatthew G. Knepley 26970034214SMatthew G. Knepley if (adjacency) { 27070034214SMatthew G. Knepley ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr); 27170034214SMatthew G. Knepley /* Get neighboring cells */ 27270034214SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 27370034214SMatthew G. Knepley PetscInt numNeighbors = PETSC_DETERMINE, n; 27470034214SMatthew G. Knepley PetscInt cellOffset = 0; 27570034214SMatthew G. Knepley 2768b0b4c70SToby Isaac ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr); 27770034214SMatthew G. Knepley /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */ 27870034214SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 27970034214SMatthew G. Knepley PetscInt cellPair[2]; 28070034214SMatthew G. Knepley PetscBool found = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE; 28170034214SMatthew G. Knepley PetscInt meetSize = 0; 28270034214SMatthew G. Knepley const PetscInt *meet = NULL; 28370034214SMatthew G. Knepley 28470034214SMatthew G. Knepley cellPair[0] = cell; cellPair[1] = neighborCells[n]; 28570034214SMatthew G. Knepley if (cellPair[0] == cellPair[1]) continue; 28670034214SMatthew G. Knepley if (!found) { 28770034214SMatthew G. Knepley ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 28870034214SMatthew G. Knepley if (meetSize) { 28970034214SMatthew G. Knepley PetscInt f; 29070034214SMatthew G. Knepley 29170034214SMatthew G. Knepley for (f = 0; f < numFaceCases; ++f) { 29270034214SMatthew G. Knepley if (numFaceVertices[f] == meetSize) { 29370034214SMatthew G. Knepley found = PETSC_TRUE; 29470034214SMatthew G. Knepley break; 29570034214SMatthew G. Knepley } 29670034214SMatthew G. Knepley } 29770034214SMatthew G. Knepley } 29870034214SMatthew G. Knepley ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 29970034214SMatthew G. Knepley } 30070034214SMatthew G. Knepley if (found) { 30170034214SMatthew G. Knepley adj[off[cell-cStart]+cellOffset] = neighborCells[n]; 30270034214SMatthew G. Knepley ++cellOffset; 30370034214SMatthew G. Knepley } 30470034214SMatthew G. Knepley } 30570034214SMatthew G. Knepley } 30670034214SMatthew G. Knepley } 30770034214SMatthew G. Knepley ierr = PetscFree(neighborCells);CHKERRQ(ierr); 30870034214SMatthew G. Knepley if (numVertices) *numVertices = numCells; 30970034214SMatthew G. Knepley if (offsets) *offsets = off; 31070034214SMatthew G. Knepley if (adjacency) *adjacency = adj; 31170034214SMatthew G. Knepley PetscFunctionReturn(0); 31270034214SMatthew G. Knepley } 31370034214SMatthew G. Knepley 31477623264SMatthew G. Knepley #undef __FUNCT__ 31577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerRegister" 31677623264SMatthew G. Knepley /*@C 31777623264SMatthew G. Knepley PetscPartitionerRegister - Adds a new PetscPartitioner implementation 31877623264SMatthew G. Knepley 31977623264SMatthew G. Knepley Not Collective 32077623264SMatthew G. Knepley 32177623264SMatthew G. Knepley Input Parameters: 32277623264SMatthew G. Knepley + name - The name of a new user-defined creation routine 32377623264SMatthew G. Knepley - create_func - The creation routine itself 32477623264SMatthew G. Knepley 32577623264SMatthew G. Knepley Notes: 32677623264SMatthew G. Knepley PetscPartitionerRegister() may be called multiple times to add several user-defined PetscPartitioners 32777623264SMatthew G. Knepley 32877623264SMatthew G. Knepley Sample usage: 32977623264SMatthew G. Knepley .vb 33077623264SMatthew G. Knepley PetscPartitionerRegister("my_part", MyPetscPartitionerCreate); 33177623264SMatthew G. Knepley .ve 33277623264SMatthew G. Knepley 33377623264SMatthew G. Knepley Then, your PetscPartitioner type can be chosen with the procedural interface via 33477623264SMatthew G. Knepley .vb 33577623264SMatthew G. Knepley PetscPartitionerCreate(MPI_Comm, PetscPartitioner *); 33677623264SMatthew G. Knepley PetscPartitionerSetType(PetscPartitioner, "my_part"); 33777623264SMatthew G. Knepley .ve 33877623264SMatthew G. Knepley or at runtime via the option 33977623264SMatthew G. Knepley .vb 34077623264SMatthew G. Knepley -petscpartitioner_type my_part 34177623264SMatthew G. Knepley .ve 34277623264SMatthew G. Knepley 34377623264SMatthew G. Knepley Level: advanced 34477623264SMatthew G. Knepley 34577623264SMatthew G. Knepley .keywords: PetscPartitioner, register 34677623264SMatthew G. Knepley .seealso: PetscPartitionerRegisterAll(), PetscPartitionerRegisterDestroy() 34777623264SMatthew G. Knepley 34877623264SMatthew G. Knepley @*/ 34977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerRegister(const char sname[], PetscErrorCode (*function)(PetscPartitioner)) 35077623264SMatthew G. Knepley { 35177623264SMatthew G. Knepley PetscErrorCode ierr; 35277623264SMatthew G. Knepley 35377623264SMatthew G. Knepley PetscFunctionBegin; 35477623264SMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscPartitionerList, sname, function);CHKERRQ(ierr); 35577623264SMatthew G. Knepley PetscFunctionReturn(0); 35677623264SMatthew G. Knepley } 35777623264SMatthew G. Knepley 35877623264SMatthew G. Knepley #undef __FUNCT__ 35977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetType" 36077623264SMatthew G. Knepley /*@C 36177623264SMatthew G. Knepley PetscPartitionerSetType - Builds a particular PetscPartitioner 36277623264SMatthew G. Knepley 36377623264SMatthew G. Knepley Collective on PetscPartitioner 36477623264SMatthew G. Knepley 36577623264SMatthew G. Knepley Input Parameters: 36677623264SMatthew G. Knepley + part - The PetscPartitioner object 36777623264SMatthew G. Knepley - name - The kind of partitioner 36877623264SMatthew G. Knepley 36977623264SMatthew G. Knepley Options Database Key: 37077623264SMatthew G. Knepley . -petscpartitioner_type <type> - Sets the PetscPartitioner type; use -help for a list of available types 37177623264SMatthew G. Knepley 37277623264SMatthew G. Knepley Level: intermediate 37377623264SMatthew G. Knepley 37477623264SMatthew G. Knepley .keywords: PetscPartitioner, set, type 37577623264SMatthew G. Knepley .seealso: PetscPartitionerGetType(), PetscPartitionerCreate() 37677623264SMatthew G. Knepley @*/ 37777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name) 37877623264SMatthew G. Knepley { 37977623264SMatthew G. Knepley PetscErrorCode (*r)(PetscPartitioner); 38077623264SMatthew G. Knepley PetscBool match; 38177623264SMatthew G. Knepley PetscErrorCode ierr; 38277623264SMatthew G. Knepley 38377623264SMatthew G. Knepley PetscFunctionBegin; 38477623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 38577623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) part, name, &match);CHKERRQ(ierr); 38677623264SMatthew G. Knepley if (match) PetscFunctionReturn(0); 38777623264SMatthew G. Knepley 3880f51fdf8SToby Isaac ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr); 38977623264SMatthew G. Knepley ierr = PetscFunctionListFind(PetscPartitionerList, name, &r);CHKERRQ(ierr); 39077623264SMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name); 39177623264SMatthew G. Knepley 39277623264SMatthew G. Knepley if (part->ops->destroy) { 39377623264SMatthew G. Knepley ierr = (*part->ops->destroy)(part);CHKERRQ(ierr); 39477623264SMatthew G. Knepley part->ops->destroy = NULL; 39577623264SMatthew G. Knepley } 39677623264SMatthew G. Knepley ierr = (*r)(part);CHKERRQ(ierr); 39777623264SMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) part, name);CHKERRQ(ierr); 39877623264SMatthew G. Knepley PetscFunctionReturn(0); 39977623264SMatthew G. Knepley } 40077623264SMatthew G. Knepley 40177623264SMatthew G. Knepley #undef __FUNCT__ 40277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerGetType" 40377623264SMatthew G. Knepley /*@C 40477623264SMatthew G. Knepley PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object. 40577623264SMatthew G. Knepley 40677623264SMatthew G. Knepley Not Collective 40777623264SMatthew G. Knepley 40877623264SMatthew G. Knepley Input Parameter: 40977623264SMatthew G. Knepley . part - The PetscPartitioner 41077623264SMatthew G. Knepley 41177623264SMatthew G. Knepley Output Parameter: 41277623264SMatthew G. Knepley . name - The PetscPartitioner type name 41377623264SMatthew G. Knepley 41477623264SMatthew G. Knepley Level: intermediate 41577623264SMatthew G. Knepley 41677623264SMatthew G. Knepley .keywords: PetscPartitioner, get, type, name 41777623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitionerCreate() 41877623264SMatthew G. Knepley @*/ 41977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name) 42077623264SMatthew G. Knepley { 42177623264SMatthew G. Knepley PetscErrorCode ierr; 42277623264SMatthew G. Knepley 42377623264SMatthew G. Knepley PetscFunctionBegin; 42477623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 42577623264SMatthew G. Knepley PetscValidPointer(name, 2); 4260f51fdf8SToby Isaac ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr); 42777623264SMatthew G. Knepley *name = ((PetscObject) part)->type_name; 42877623264SMatthew G. Knepley PetscFunctionReturn(0); 42977623264SMatthew G. Knepley } 43077623264SMatthew G. Knepley 43177623264SMatthew G. Knepley #undef __FUNCT__ 43277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView" 43377623264SMatthew G. Knepley /*@C 43477623264SMatthew G. Knepley PetscPartitionerView - Views a PetscPartitioner 43577623264SMatthew G. Knepley 43677623264SMatthew G. Knepley Collective on PetscPartitioner 43777623264SMatthew G. Knepley 43877623264SMatthew G. Knepley Input Parameter: 43977623264SMatthew G. Knepley + part - the PetscPartitioner object to view 44077623264SMatthew G. Knepley - v - the viewer 44177623264SMatthew G. Knepley 44277623264SMatthew G. Knepley Level: developer 44377623264SMatthew G. Knepley 44477623264SMatthew G. Knepley .seealso: PetscPartitionerDestroy() 44577623264SMatthew G. Knepley @*/ 44677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v) 44777623264SMatthew G. Knepley { 44877623264SMatthew G. Knepley PetscErrorCode ierr; 44977623264SMatthew G. Knepley 45077623264SMatthew G. Knepley PetscFunctionBegin; 45177623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 45277623264SMatthew G. Knepley if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) part), &v);CHKERRQ(ierr);} 45377623264SMatthew G. Knepley if (part->ops->view) {ierr = (*part->ops->view)(part, v);CHKERRQ(ierr);} 45477623264SMatthew G. Knepley PetscFunctionReturn(0); 45577623264SMatthew G. Knepley } 45677623264SMatthew G. Knepley 45777623264SMatthew G. Knepley #undef __FUNCT__ 45877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetTypeFromOptions_Internal" 45977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner part) 46077623264SMatthew G. Knepley { 46177623264SMatthew G. Knepley const char *defaultType; 46277623264SMatthew G. Knepley char name[256]; 46377623264SMatthew G. Knepley PetscBool flg; 46477623264SMatthew G. Knepley PetscErrorCode ierr; 46577623264SMatthew G. Knepley 46677623264SMatthew G. Knepley PetscFunctionBegin; 46777623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 46877623264SMatthew G. Knepley if (!((PetscObject) part)->type_name) defaultType = PETSCPARTITIONERCHACO; 46977623264SMatthew G. Knepley else defaultType = ((PetscObject) part)->type_name; 4700f51fdf8SToby Isaac ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr); 47177623264SMatthew G. Knepley 47277623264SMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr); 47377623264SMatthew G. Knepley ierr = PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, defaultType, name, 256, &flg);CHKERRQ(ierr); 47477623264SMatthew G. Knepley if (flg) { 47577623264SMatthew G. Knepley ierr = PetscPartitionerSetType(part, name);CHKERRQ(ierr); 47677623264SMatthew G. Knepley } else if (!((PetscObject) part)->type_name) { 47777623264SMatthew G. Knepley ierr = PetscPartitionerSetType(part, defaultType);CHKERRQ(ierr); 47877623264SMatthew G. Knepley } 47977623264SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 48077623264SMatthew G. Knepley PetscFunctionReturn(0); 48177623264SMatthew G. Knepley } 48277623264SMatthew G. Knepley 48377623264SMatthew G. Knepley #undef __FUNCT__ 48477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetFromOptions" 48577623264SMatthew G. Knepley /*@ 48677623264SMatthew G. Knepley PetscPartitionerSetFromOptions - sets parameters in a PetscPartitioner from the options database 48777623264SMatthew G. Knepley 48877623264SMatthew G. Knepley Collective on PetscPartitioner 48977623264SMatthew G. Knepley 49077623264SMatthew G. Knepley Input Parameter: 49177623264SMatthew G. Knepley . part - the PetscPartitioner object to set options for 49277623264SMatthew G. Knepley 49377623264SMatthew G. Knepley Level: developer 49477623264SMatthew G. Knepley 49577623264SMatthew G. Knepley .seealso: PetscPartitionerView() 49677623264SMatthew G. Knepley @*/ 49777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part) 49877623264SMatthew G. Knepley { 49977623264SMatthew G. Knepley PetscErrorCode ierr; 50077623264SMatthew G. Knepley 50177623264SMatthew G. Knepley PetscFunctionBegin; 50277623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 50377623264SMatthew G. Knepley ierr = PetscPartitionerSetTypeFromOptions_Internal(part);CHKERRQ(ierr); 50477623264SMatthew G. Knepley 50577623264SMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr); 50677623264SMatthew G. Knepley if (part->ops->setfromoptions) {ierr = (*part->ops->setfromoptions)(part);CHKERRQ(ierr);} 50777623264SMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 508*0633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) part);CHKERRQ(ierr); 50977623264SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 51077623264SMatthew G. Knepley ierr = PetscPartitionerViewFromOptions(part, NULL, "-petscpartitioner_view");CHKERRQ(ierr); 51177623264SMatthew G. Knepley PetscFunctionReturn(0); 51277623264SMatthew G. Knepley } 51377623264SMatthew G. Knepley 51477623264SMatthew G. Knepley #undef __FUNCT__ 51577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetUp" 51677623264SMatthew G. Knepley /*@C 51777623264SMatthew G. Knepley PetscPartitionerSetUp - Construct data structures for the PetscPartitioner 51877623264SMatthew G. Knepley 51977623264SMatthew G. Knepley Collective on PetscPartitioner 52077623264SMatthew G. Knepley 52177623264SMatthew G. Knepley Input Parameter: 52277623264SMatthew G. Knepley . part - the PetscPartitioner object to setup 52377623264SMatthew G. Knepley 52477623264SMatthew G. Knepley Level: developer 52577623264SMatthew G. Knepley 52677623264SMatthew G. Knepley .seealso: PetscPartitionerView(), PetscPartitionerDestroy() 52777623264SMatthew G. Knepley @*/ 52877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part) 52977623264SMatthew G. Knepley { 53077623264SMatthew G. Knepley PetscErrorCode ierr; 53177623264SMatthew G. Knepley 53277623264SMatthew G. Knepley PetscFunctionBegin; 53377623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 53477623264SMatthew G. Knepley if (part->ops->setup) {ierr = (*part->ops->setup)(part);CHKERRQ(ierr);} 53577623264SMatthew G. Knepley PetscFunctionReturn(0); 53677623264SMatthew G. Knepley } 53777623264SMatthew G. Knepley 53877623264SMatthew G. Knepley #undef __FUNCT__ 53977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy" 54077623264SMatthew G. Knepley /*@ 54177623264SMatthew G. Knepley PetscPartitionerDestroy - Destroys a PetscPartitioner object 54277623264SMatthew G. Knepley 54377623264SMatthew G. Knepley Collective on PetscPartitioner 54477623264SMatthew G. Knepley 54577623264SMatthew G. Knepley Input Parameter: 54677623264SMatthew G. Knepley . part - the PetscPartitioner object to destroy 54777623264SMatthew G. Knepley 54877623264SMatthew G. Knepley Level: developer 54977623264SMatthew G. Knepley 55077623264SMatthew G. Knepley .seealso: PetscPartitionerView() 55177623264SMatthew G. Knepley @*/ 55277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part) 55377623264SMatthew G. Knepley { 55477623264SMatthew G. Knepley PetscErrorCode ierr; 55577623264SMatthew G. Knepley 55677623264SMatthew G. Knepley PetscFunctionBegin; 55777623264SMatthew G. Knepley if (!*part) PetscFunctionReturn(0); 55877623264SMatthew G. Knepley PetscValidHeaderSpecific((*part), PETSCPARTITIONER_CLASSID, 1); 55977623264SMatthew G. Knepley 56077623264SMatthew G. Knepley if (--((PetscObject)(*part))->refct > 0) {*part = 0; PetscFunctionReturn(0);} 56177623264SMatthew G. Knepley ((PetscObject) (*part))->refct = 0; 56277623264SMatthew G. Knepley 56377623264SMatthew G. Knepley if ((*part)->ops->destroy) {ierr = (*(*part)->ops->destroy)(*part);CHKERRQ(ierr);} 56477623264SMatthew G. Knepley ierr = PetscHeaderDestroy(part);CHKERRQ(ierr); 56577623264SMatthew G. Knepley PetscFunctionReturn(0); 56677623264SMatthew G. Knepley } 56777623264SMatthew G. Knepley 56877623264SMatthew G. Knepley #undef __FUNCT__ 56977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate" 57077623264SMatthew G. Knepley /*@ 57177623264SMatthew G. Knepley PetscPartitionerCreate - Creates an empty PetscPartitioner object. The type can then be set with PetscPartitionerSetType(). 57277623264SMatthew G. Knepley 57377623264SMatthew G. Knepley Collective on MPI_Comm 57477623264SMatthew G. Knepley 57577623264SMatthew G. Knepley Input Parameter: 57677623264SMatthew G. Knepley . comm - The communicator for the PetscPartitioner object 57777623264SMatthew G. Knepley 57877623264SMatthew G. Knepley Output Parameter: 57977623264SMatthew G. Knepley . part - The PetscPartitioner object 58077623264SMatthew G. Knepley 58177623264SMatthew G. Knepley Level: beginner 58277623264SMatthew G. Knepley 583555a9cf8SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PETSCPARTITIONERCHACO, PETSCPARTITIONERPARMETIS, PETSCPARTITIONERSHELL, PETSCPARTITIONERSIMPLE 58477623264SMatthew G. Knepley @*/ 58577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part) 58677623264SMatthew G. Knepley { 58777623264SMatthew G. Knepley PetscPartitioner p; 58877623264SMatthew G. Knepley PetscErrorCode ierr; 58977623264SMatthew G. Knepley 59077623264SMatthew G. Knepley PetscFunctionBegin; 59177623264SMatthew G. Knepley PetscValidPointer(part, 2); 59277623264SMatthew G. Knepley *part = NULL; 59377623264SMatthew G. Knepley ierr = PetscFVInitializePackage();CHKERRQ(ierr); 59477623264SMatthew G. Knepley 59573107ff1SLisandro Dalcin ierr = PetscHeaderCreate(p, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView);CHKERRQ(ierr); 59677623264SMatthew G. Knepley 59777623264SMatthew G. Knepley *part = p; 59877623264SMatthew G. Knepley PetscFunctionReturn(0); 59977623264SMatthew G. Knepley } 60077623264SMatthew G. Knepley 60177623264SMatthew G. Knepley #undef __FUNCT__ 60277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition" 60377623264SMatthew G. Knepley /*@ 60477623264SMatthew G. Knepley PetscPartitionerPartition - Create a non-overlapping partition of the cells in the mesh 60577623264SMatthew G. Knepley 60677623264SMatthew G. Knepley Collective on DM 60777623264SMatthew G. Knepley 60877623264SMatthew G. Knepley Input Parameters: 60977623264SMatthew G. Knepley + part - The PetscPartitioner 610f8987ae8SMichael Lange - dm - The mesh DM 61177623264SMatthew G. Knepley 61277623264SMatthew G. Knepley Output Parameters: 61377623264SMatthew G. Knepley + partSection - The PetscSection giving the division of points by partition 614f8987ae8SMichael Lange - partition - The list of points by partition 61577623264SMatthew G. Knepley 61677623264SMatthew G. Knepley Note: Instead of cells, points at a given height can be partitioned by calling PetscPartitionerSetPointHeight() 61777623264SMatthew G. Knepley 61877623264SMatthew G. Knepley Level: developer 61977623264SMatthew G. Knepley 62077623264SMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerSetPointHeight(), PetscPartitionerCreate() 6214b15ede2SMatthew G. Knepley @*/ 622f8987ae8SMichael Lange PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, DM dm, PetscSection partSection, IS *partition) 62377623264SMatthew G. Knepley { 62477623264SMatthew G. Knepley PetscMPIInt size; 62577623264SMatthew G. Knepley PetscErrorCode ierr; 62677623264SMatthew G. Knepley 62777623264SMatthew G. Knepley PetscFunctionBegin; 62877623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 62977623264SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 63077623264SMatthew G. Knepley PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 4); 63177623264SMatthew G. Knepley PetscValidPointer(partition, 5); 63277623264SMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) part), &size);CHKERRQ(ierr); 63377623264SMatthew G. Knepley if (size == 1) { 63477623264SMatthew G. Knepley PetscInt *points; 63577623264SMatthew G. Knepley PetscInt cStart, cEnd, c; 63677623264SMatthew G. Knepley 63777623264SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr); 63877623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, size);CHKERRQ(ierr); 63977623264SMatthew G. Knepley ierr = PetscSectionSetDof(partSection, 0, cEnd-cStart);CHKERRQ(ierr); 64077623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 64177623264SMatthew G. Knepley ierr = PetscMalloc1(cEnd-cStart, &points);CHKERRQ(ierr); 64277623264SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) points[c] = c; 64377623264SMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 64477623264SMatthew G. Knepley PetscFunctionReturn(0); 64577623264SMatthew G. Knepley } 64677623264SMatthew G. Knepley if (part->height == 0) { 64777623264SMatthew G. Knepley PetscInt numVertices; 64877623264SMatthew G. Knepley PetscInt *start = NULL; 64977623264SMatthew G. Knepley PetscInt *adjacency = NULL; 65077623264SMatthew G. Knepley 651532c4e7dSMichael Lange ierr = DMPlexCreatePartitionerGraph(dm, 0, &numVertices, &start, &adjacency);CHKERRQ(ierr); 65277623264SMatthew G. Knepley if (!part->ops->partition) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONGSTATE, "PetscPartitioner has no type"); 65377623264SMatthew G. Knepley ierr = (*part->ops->partition)(part, dm, size, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr); 65477623264SMatthew G. Knepley ierr = PetscFree(start);CHKERRQ(ierr); 65577623264SMatthew G. Knepley ierr = PetscFree(adjacency);CHKERRQ(ierr); 65677623264SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_OUTOFRANGE, "Invalid height %D for points to partition", part->height); 65777623264SMatthew G. Knepley PetscFunctionReturn(0); 65877623264SMatthew G. Knepley } 65977623264SMatthew G. Knepley 66077623264SMatthew G. Knepley #undef __FUNCT__ 66177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Shell" 66277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part) 66377623264SMatthew G. Knepley { 66477623264SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 66577623264SMatthew G. Knepley PetscErrorCode ierr; 66677623264SMatthew G. Knepley 66777623264SMatthew G. Knepley PetscFunctionBegin; 66877623264SMatthew G. Knepley ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr); 66977623264SMatthew G. Knepley ierr = ISDestroy(&p->partition);CHKERRQ(ierr); 67077623264SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 67177623264SMatthew G. Knepley PetscFunctionReturn(0); 67277623264SMatthew G. Knepley } 67377623264SMatthew G. Knepley 67477623264SMatthew G. Knepley #undef __FUNCT__ 67577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Shell_Ascii" 67677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell_Ascii(PetscPartitioner part, PetscViewer viewer) 67777623264SMatthew G. Knepley { 67877623264SMatthew G. Knepley PetscViewerFormat format; 67977623264SMatthew G. Knepley PetscErrorCode ierr; 68077623264SMatthew G. Knepley 68177623264SMatthew G. Knepley PetscFunctionBegin; 68277623264SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 68377623264SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Shell Graph Partitioner:\n");CHKERRQ(ierr); 68477623264SMatthew G. Knepley PetscFunctionReturn(0); 68577623264SMatthew G. Knepley } 68677623264SMatthew G. Knepley 68777623264SMatthew G. Knepley #undef __FUNCT__ 68877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Shell" 68977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer) 69077623264SMatthew G. Knepley { 69177623264SMatthew G. Knepley PetscBool iascii; 69277623264SMatthew G. Knepley PetscErrorCode ierr; 69377623264SMatthew G. Knepley 69477623264SMatthew G. Knepley PetscFunctionBegin; 69577623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 69677623264SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 69777623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 69877623264SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_Shell_Ascii(part, viewer);CHKERRQ(ierr);} 69977623264SMatthew G. Knepley PetscFunctionReturn(0); 70077623264SMatthew G. Knepley } 70177623264SMatthew G. Knepley 70277623264SMatthew G. Knepley #undef __FUNCT__ 70377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Shell" 70477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 70577623264SMatthew G. Knepley { 70677623264SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 70777623264SMatthew G. Knepley PetscInt np; 70877623264SMatthew G. Knepley PetscErrorCode ierr; 70977623264SMatthew G. Knepley 71077623264SMatthew G. Knepley PetscFunctionBegin; 71177623264SMatthew G. Knepley ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr); 71277623264SMatthew G. Knepley if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np); 71377623264SMatthew G. Knepley ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr); 71477623264SMatthew G. Knepley if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np); 7155680f57bSMatthew G. Knepley ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr); 71677623264SMatthew G. Knepley *partition = p->partition; 71777623264SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr); 71877623264SMatthew G. Knepley PetscFunctionReturn(0); 71977623264SMatthew G. Knepley } 72077623264SMatthew G. Knepley 72177623264SMatthew G. Knepley #undef __FUNCT__ 72277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Shell" 72377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part) 72477623264SMatthew G. Knepley { 72577623264SMatthew G. Knepley PetscFunctionBegin; 72677623264SMatthew G. Knepley part->ops->view = PetscPartitionerView_Shell; 72777623264SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_Shell; 72877623264SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_Shell; 72977623264SMatthew G. Knepley PetscFunctionReturn(0); 73077623264SMatthew G. Knepley } 73177623264SMatthew G. Knepley 73277623264SMatthew G. Knepley /*MC 73377623264SMatthew G. Knepley PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object 73477623264SMatthew G. Knepley 73577623264SMatthew G. Knepley Level: intermediate 73677623264SMatthew G. Knepley 73777623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 73877623264SMatthew G. Knepley M*/ 73977623264SMatthew G. Knepley 74077623264SMatthew G. Knepley #undef __FUNCT__ 74177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Shell" 74277623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part) 74377623264SMatthew G. Knepley { 74477623264SMatthew G. Knepley PetscPartitioner_Shell *p; 74577623264SMatthew G. Knepley PetscErrorCode ierr; 74677623264SMatthew G. Knepley 74777623264SMatthew G. Knepley PetscFunctionBegin; 74877623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 74977623264SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 75077623264SMatthew G. Knepley part->data = p; 75177623264SMatthew G. Knepley 75277623264SMatthew G. Knepley ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr); 75377623264SMatthew G. Knepley PetscFunctionReturn(0); 75477623264SMatthew G. Knepley } 75577623264SMatthew G. Knepley 75677623264SMatthew G. Knepley #undef __FUNCT__ 7575680f57bSMatthew G. Knepley #define __FUNCT__ "PetscPartitionerShellSetPartition" 7585680f57bSMatthew G. Knepley /*@C 7595680f57bSMatthew G. Knepley PetscPartitionerShellSetPartition - Set an artifical partition for a mesh 7605680f57bSMatthew G. Knepley 761b7e49471SLawrence Mitchell Collective on PART 7625680f57bSMatthew G. Knepley 7635680f57bSMatthew G. Knepley Input Parameters: 7645680f57bSMatthew G. Knepley + part - The PetscPartitioner 765b7e49471SLawrence Mitchell . numProcs - The number of partitions 766b7e49471SLawrence Mitchell . sizes - array of size numProcs (or NULL) providing the number of points in each partition 767b7e49471SLawrence Mitchell - points - array of size sum(sizes) (may be NULL iff sizes is NULL) providing the partition each point belongs to 7685680f57bSMatthew G. Knepley 7695680f57bSMatthew G. Knepley Level: developer 7705680f57bSMatthew G. Knepley 771b7e49471SLawrence Mitchell Notes: 772b7e49471SLawrence Mitchell 773b7e49471SLawrence Mitchell It is safe to free the sizes and points arrays after use in this routine. 774b7e49471SLawrence Mitchell 7755680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate() 7765680f57bSMatthew G. Knepley @*/ 7775680f57bSMatthew G. Knepley PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt numProcs, const PetscInt sizes[], const PetscInt points[]) 7785680f57bSMatthew G. Knepley { 7795680f57bSMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 7805680f57bSMatthew G. Knepley PetscInt proc, numPoints; 7815680f57bSMatthew G. Knepley PetscErrorCode ierr; 7825680f57bSMatthew G. Knepley 7835680f57bSMatthew G. Knepley PetscFunctionBegin; 7845680f57bSMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 7855680f57bSMatthew G. Knepley if (sizes) {PetscValidPointer(sizes, 3);} 7865680f57bSMatthew G. Knepley if (sizes) {PetscValidPointer(points, 4);} 7875680f57bSMatthew G. Knepley ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr); 7885680f57bSMatthew G. Knepley ierr = ISDestroy(&p->partition);CHKERRQ(ierr); 7895680f57bSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr); 7905680f57bSMatthew G. Knepley ierr = PetscSectionSetChart(p->section, 0, numProcs);CHKERRQ(ierr); 7915680f57bSMatthew G. Knepley if (sizes) { 7925680f57bSMatthew G. Knepley for (proc = 0; proc < numProcs; ++proc) { 7935680f57bSMatthew G. Knepley ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr); 7945680f57bSMatthew G. Knepley } 7955680f57bSMatthew G. Knepley } 7965680f57bSMatthew G. Knepley ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr); 7975680f57bSMatthew G. Knepley ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr); 7985680f57bSMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr); 7995680f57bSMatthew G. Knepley PetscFunctionReturn(0); 8005680f57bSMatthew G. Knepley } 8015680f57bSMatthew G. Knepley 8025680f57bSMatthew G. Knepley #undef __FUNCT__ 803555a9cf8SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Simple" 804555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part) 805555a9cf8SMatthew G. Knepley { 806555a9cf8SMatthew G. Knepley PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 807555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 808555a9cf8SMatthew G. Knepley 809555a9cf8SMatthew G. Knepley PetscFunctionBegin; 810555a9cf8SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 811555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 812555a9cf8SMatthew G. Knepley } 813555a9cf8SMatthew G. Knepley 814555a9cf8SMatthew G. Knepley #undef __FUNCT__ 815555a9cf8SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Simple_Ascii" 816555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Simple_Ascii(PetscPartitioner part, PetscViewer viewer) 817555a9cf8SMatthew G. Knepley { 818555a9cf8SMatthew G. Knepley PetscViewerFormat format; 819555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 820555a9cf8SMatthew G. Knepley 821555a9cf8SMatthew G. Knepley PetscFunctionBegin; 822555a9cf8SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 823555a9cf8SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Simple Graph Partitioner:\n");CHKERRQ(ierr); 824555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 825555a9cf8SMatthew G. Knepley } 826555a9cf8SMatthew G. Knepley 827555a9cf8SMatthew G. Knepley #undef __FUNCT__ 828555a9cf8SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Simple" 829555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer) 830555a9cf8SMatthew G. Knepley { 831555a9cf8SMatthew G. Knepley PetscBool iascii; 832555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 833555a9cf8SMatthew G. Knepley 834555a9cf8SMatthew G. Knepley PetscFunctionBegin; 835555a9cf8SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 836555a9cf8SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 837555a9cf8SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 838555a9cf8SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_Simple_Ascii(part, viewer);CHKERRQ(ierr);} 839555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 840555a9cf8SMatthew G. Knepley } 841555a9cf8SMatthew G. Knepley 842555a9cf8SMatthew G. Knepley #undef __FUNCT__ 843555a9cf8SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Simple" 844555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 845555a9cf8SMatthew G. Knepley { 846555a9cf8SMatthew G. Knepley PetscInt np; 847555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 848555a9cf8SMatthew G. Knepley 849555a9cf8SMatthew G. Knepley PetscFunctionBegin; 850555a9cf8SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 851555a9cf8SMatthew G. Knepley for (np = 0; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np));CHKERRQ(ierr);} 852555a9cf8SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 853555a9cf8SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr); 854555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 855555a9cf8SMatthew G. Knepley } 856555a9cf8SMatthew G. Knepley 857555a9cf8SMatthew G. Knepley #undef __FUNCT__ 858555a9cf8SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Simple" 859555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part) 860555a9cf8SMatthew G. Knepley { 861555a9cf8SMatthew G. Knepley PetscFunctionBegin; 862555a9cf8SMatthew G. Knepley part->ops->view = PetscPartitionerView_Simple; 863555a9cf8SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_Simple; 864555a9cf8SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_Simple; 865555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 866555a9cf8SMatthew G. Knepley } 867555a9cf8SMatthew G. Knepley 868555a9cf8SMatthew G. Knepley /*MC 869555a9cf8SMatthew G. Knepley PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object 870555a9cf8SMatthew G. Knepley 871555a9cf8SMatthew G. Knepley Level: intermediate 872555a9cf8SMatthew G. Knepley 873555a9cf8SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 874555a9cf8SMatthew G. Knepley M*/ 875555a9cf8SMatthew G. Knepley 876555a9cf8SMatthew G. Knepley #undef __FUNCT__ 877555a9cf8SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Simple" 878555a9cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part) 879555a9cf8SMatthew G. Knepley { 880555a9cf8SMatthew G. Knepley PetscPartitioner_Simple *p; 881555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 882555a9cf8SMatthew G. Knepley 883555a9cf8SMatthew G. Knepley PetscFunctionBegin; 884555a9cf8SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 885555a9cf8SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 886555a9cf8SMatthew G. Knepley part->data = p; 887555a9cf8SMatthew G. Knepley 888555a9cf8SMatthew G. Knepley ierr = PetscPartitionerInitialize_Simple(part);CHKERRQ(ierr); 889555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 890555a9cf8SMatthew G. Knepley } 891555a9cf8SMatthew G. Knepley 892555a9cf8SMatthew G. Knepley #undef __FUNCT__ 89377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Chaco" 89477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part) 89577623264SMatthew G. Knepley { 89677623264SMatthew G. Knepley PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data; 89777623264SMatthew G. Knepley PetscErrorCode ierr; 89877623264SMatthew G. Knepley 89977623264SMatthew G. Knepley PetscFunctionBegin; 90077623264SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 90177623264SMatthew G. Knepley PetscFunctionReturn(0); 90277623264SMatthew G. Knepley } 90377623264SMatthew G. Knepley 90477623264SMatthew G. Knepley #undef __FUNCT__ 90577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Chaco_Ascii" 90677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco_Ascii(PetscPartitioner part, PetscViewer viewer) 90777623264SMatthew G. Knepley { 90877623264SMatthew G. Knepley PetscViewerFormat format; 90977623264SMatthew G. Knepley PetscErrorCode ierr; 91077623264SMatthew G. Knepley 91177623264SMatthew G. Knepley PetscFunctionBegin; 91277623264SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 91377623264SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Chaco Graph Partitioner:\n");CHKERRQ(ierr); 91477623264SMatthew G. Knepley PetscFunctionReturn(0); 91577623264SMatthew G. Knepley } 91677623264SMatthew G. Knepley 91777623264SMatthew G. Knepley #undef __FUNCT__ 91877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Chaco" 91977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer) 92077623264SMatthew G. Knepley { 92177623264SMatthew G. Knepley PetscBool iascii; 92277623264SMatthew G. Knepley PetscErrorCode ierr; 92377623264SMatthew G. Knepley 92477623264SMatthew G. Knepley PetscFunctionBegin; 92577623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 92677623264SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 92777623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 92877623264SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_Chaco_Ascii(part, viewer);CHKERRQ(ierr);} 92977623264SMatthew G. Knepley PetscFunctionReturn(0); 93077623264SMatthew G. Knepley } 93177623264SMatthew G. Knepley 93270034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO) 93370034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H) 93470034214SMatthew G. Knepley #include <unistd.h> 93570034214SMatthew G. Knepley #endif 93670034214SMatthew G. Knepley /* Chaco does not have an include file */ 93770034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts, 93870034214SMatthew G. Knepley float *ewgts, float *x, float *y, float *z, char *outassignname, 93970034214SMatthew G. Knepley char *outfilename, short *assignment, int architecture, int ndims_tot, 94070034214SMatthew G. Knepley int mesh_dims[3], double *goal, int global_method, int local_method, 94170034214SMatthew G. Knepley int rqi_flag, int vmax, int ndims, double eigtol, long seed); 94270034214SMatthew G. Knepley 94370034214SMatthew G. Knepley extern int FREE_GRAPH; 94477623264SMatthew G. Knepley #endif 94570034214SMatthew G. Knepley 94670034214SMatthew G. Knepley #undef __FUNCT__ 94777623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Chaco" 94877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 94970034214SMatthew G. Knepley { 95077623264SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO) 95170034214SMatthew G. Knepley enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3}; 95270034214SMatthew G. Knepley MPI_Comm comm; 95370034214SMatthew G. Knepley int nvtxs = numVertices; /* number of vertices in full graph */ 95470034214SMatthew G. Knepley int *vwgts = NULL; /* weights for all vertices */ 95570034214SMatthew G. Knepley float *ewgts = NULL; /* weights for all edges */ 95670034214SMatthew G. Knepley float *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */ 95770034214SMatthew G. Knepley char *outassignname = NULL; /* name of assignment output file */ 95870034214SMatthew G. Knepley char *outfilename = NULL; /* output file name */ 95970034214SMatthew G. Knepley int architecture = 1; /* 0 => hypercube, d => d-dimensional mesh */ 96070034214SMatthew G. Knepley int ndims_tot = 0; /* total number of cube dimensions to divide */ 96170034214SMatthew G. Knepley int mesh_dims[3]; /* dimensions of mesh of processors */ 96270034214SMatthew G. Knepley double *goal = NULL; /* desired set sizes for each set */ 96370034214SMatthew G. Knepley int global_method = 1; /* global partitioning algorithm */ 96470034214SMatthew G. Knepley int local_method = 1; /* local partitioning algorithm */ 96570034214SMatthew G. Knepley int rqi_flag = 0; /* should I use RQI/Symmlq eigensolver? */ 96670034214SMatthew G. Knepley int vmax = 200; /* how many vertices to coarsen down to? */ 96770034214SMatthew G. Knepley int ndims = 1; /* number of eigenvectors (2^d sets) */ 96870034214SMatthew G. Knepley double eigtol = 0.001; /* tolerance on eigenvectors */ 96970034214SMatthew G. Knepley long seed = 123636512; /* for random graph mutations */ 97070034214SMatthew G. Knepley short int *assignment; /* Output partition */ 97170034214SMatthew G. Knepley int fd_stdout, fd_pipe[2]; 97270034214SMatthew G. Knepley PetscInt *points; 97370034214SMatthew G. Knepley int i, v, p; 97470034214SMatthew G. Knepley PetscErrorCode ierr; 97570034214SMatthew G. Knepley 97670034214SMatthew G. Knepley PetscFunctionBegin; 97770034214SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 97870034214SMatthew G. Knepley if (!numVertices) { 97977623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 98077623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 98170034214SMatthew G. Knepley ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 98270034214SMatthew G. Knepley PetscFunctionReturn(0); 98370034214SMatthew G. Knepley } 98470034214SMatthew G. Knepley FREE_GRAPH = 0; /* Do not let Chaco free my memory */ 98570034214SMatthew G. Knepley for (i = 0; i < start[numVertices]; ++i) ++adjacency[i]; 98670034214SMatthew G. Knepley 98770034214SMatthew G. Knepley if (global_method == INERTIAL_METHOD) { 98870034214SMatthew G. Knepley /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */ 98970034214SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported"); 99070034214SMatthew G. Knepley } 99177623264SMatthew G. Knepley mesh_dims[0] = nparts; 99270034214SMatthew G. Knepley mesh_dims[1] = 1; 99370034214SMatthew G. Knepley mesh_dims[2] = 1; 99470034214SMatthew G. Knepley ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr); 99570034214SMatthew G. Knepley /* Chaco outputs to stdout. We redirect this to a buffer. */ 99670034214SMatthew G. Knepley /* TODO: check error codes for UNIX calls */ 99770034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H) 99870034214SMatthew G. Knepley { 99970034214SMatthew G. Knepley int piperet; 100070034214SMatthew G. Knepley piperet = pipe(fd_pipe); 100170034214SMatthew G. Knepley if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe"); 100270034214SMatthew G. Knepley fd_stdout = dup(1); 100370034214SMatthew G. Knepley close(1); 100470034214SMatthew G. Knepley dup2(fd_pipe[1], 1); 100570034214SMatthew G. Knepley } 100670034214SMatthew G. Knepley #endif 100770034214SMatthew G. Knepley ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename, 100870034214SMatthew G. Knepley assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag, 100970034214SMatthew G. Knepley vmax, ndims, eigtol, seed); 101070034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H) 101170034214SMatthew G. Knepley { 101270034214SMatthew G. Knepley char msgLog[10000]; 101370034214SMatthew G. Knepley int count; 101470034214SMatthew G. Knepley 101570034214SMatthew G. Knepley fflush(stdout); 101670034214SMatthew G. Knepley count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char)); 101770034214SMatthew G. Knepley if (count < 0) count = 0; 101870034214SMatthew G. Knepley msgLog[count] = 0; 101970034214SMatthew G. Knepley close(1); 102070034214SMatthew G. Knepley dup2(fd_stdout, 1); 102170034214SMatthew G. Knepley close(fd_stdout); 102270034214SMatthew G. Knepley close(fd_pipe[0]); 102370034214SMatthew G. Knepley close(fd_pipe[1]); 102470034214SMatthew G. Knepley if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog); 102570034214SMatthew G. Knepley } 102670034214SMatthew G. Knepley #endif 102770034214SMatthew G. Knepley /* Convert to PetscSection+IS */ 102877623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 102970034214SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) { 103077623264SMatthew G. Knepley ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr); 103170034214SMatthew G. Knepley } 103277623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 103370034214SMatthew G. Knepley ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr); 103477623264SMatthew G. Knepley for (p = 0, i = 0; p < nparts; ++p) { 103570034214SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) { 103670034214SMatthew G. Knepley if (assignment[v] == p) points[i++] = v; 103770034214SMatthew G. Knepley } 103870034214SMatthew G. Knepley } 103970034214SMatthew G. Knepley if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs); 104070034214SMatthew G. Knepley ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 104170034214SMatthew G. Knepley if (global_method == INERTIAL_METHOD) { 104270034214SMatthew G. Knepley /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */ 104370034214SMatthew G. Knepley } 104470034214SMatthew G. Knepley ierr = PetscFree(assignment);CHKERRQ(ierr); 104570034214SMatthew G. Knepley for (i = 0; i < start[numVertices]; ++i) --adjacency[i]; 104670034214SMatthew G. Knepley PetscFunctionReturn(0); 104777623264SMatthew G. Knepley #else 104877623264SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco."); 104970034214SMatthew G. Knepley #endif 105077623264SMatthew G. Knepley } 105177623264SMatthew G. Knepley 105277623264SMatthew G. Knepley #undef __FUNCT__ 105377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Chaco" 105477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part) 105577623264SMatthew G. Knepley { 105677623264SMatthew G. Knepley PetscFunctionBegin; 105777623264SMatthew G. Knepley part->ops->view = PetscPartitionerView_Chaco; 105877623264SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_Chaco; 105977623264SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_Chaco; 106077623264SMatthew G. Knepley PetscFunctionReturn(0); 106177623264SMatthew G. Knepley } 106277623264SMatthew G. Knepley 106377623264SMatthew G. Knepley /*MC 106477623264SMatthew G. Knepley PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library 106577623264SMatthew G. Knepley 106677623264SMatthew G. Knepley Level: intermediate 106777623264SMatthew G. Knepley 106877623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 106977623264SMatthew G. Knepley M*/ 107077623264SMatthew G. Knepley 107177623264SMatthew G. Knepley #undef __FUNCT__ 107277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Chaco" 107377623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part) 107477623264SMatthew G. Knepley { 107577623264SMatthew G. Knepley PetscPartitioner_Chaco *p; 107677623264SMatthew G. Knepley PetscErrorCode ierr; 107777623264SMatthew G. Knepley 107877623264SMatthew G. Knepley PetscFunctionBegin; 107977623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 108077623264SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 108177623264SMatthew G. Knepley part->data = p; 108277623264SMatthew G. Knepley 108377623264SMatthew G. Knepley ierr = PetscPartitionerInitialize_Chaco(part);CHKERRQ(ierr); 108477623264SMatthew G. Knepley ierr = PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionercite);CHKERRQ(ierr); 108577623264SMatthew G. Knepley PetscFunctionReturn(0); 108677623264SMatthew G. Knepley } 108777623264SMatthew G. Knepley 108877623264SMatthew G. Knepley #undef __FUNCT__ 108977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_ParMetis" 109077623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part) 109177623264SMatthew G. Knepley { 109277623264SMatthew G. Knepley PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *) part->data; 109377623264SMatthew G. Knepley PetscErrorCode ierr; 109477623264SMatthew G. Knepley 109577623264SMatthew G. Knepley PetscFunctionBegin; 109677623264SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 109777623264SMatthew G. Knepley PetscFunctionReturn(0); 109877623264SMatthew G. Knepley } 109977623264SMatthew G. Knepley 110077623264SMatthew G. Knepley #undef __FUNCT__ 110177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_ParMetis_Ascii" 110277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis_Ascii(PetscPartitioner part, PetscViewer viewer) 110377623264SMatthew G. Knepley { 110477623264SMatthew G. Knepley PetscViewerFormat format; 110577623264SMatthew G. Knepley PetscErrorCode ierr; 110677623264SMatthew G. Knepley 110777623264SMatthew G. Knepley PetscFunctionBegin; 110877623264SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 110977623264SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "ParMetis Graph Partitioner:\n");CHKERRQ(ierr); 111077623264SMatthew G. Knepley PetscFunctionReturn(0); 111177623264SMatthew G. Knepley } 111277623264SMatthew G. Knepley 111377623264SMatthew G. Knepley #undef __FUNCT__ 111477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_ParMetis" 111577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer) 111677623264SMatthew G. Knepley { 111777623264SMatthew G. Knepley PetscBool iascii; 111877623264SMatthew G. Knepley PetscErrorCode ierr; 111977623264SMatthew G. Knepley 112077623264SMatthew G. Knepley PetscFunctionBegin; 112177623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 112277623264SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 112377623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 112477623264SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_ParMetis_Ascii(part, viewer);CHKERRQ(ierr);} 112577623264SMatthew G. Knepley PetscFunctionReturn(0); 112677623264SMatthew G. Knepley } 112770034214SMatthew G. Knepley 112870034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS) 112970034214SMatthew G. Knepley #include <parmetis.h> 113077623264SMatthew G. Knepley #endif 113170034214SMatthew G. Knepley 113270034214SMatthew G. Knepley #undef __FUNCT__ 113377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_ParMetis" 113477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 113570034214SMatthew G. Knepley { 113677623264SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS) 113770034214SMatthew G. Knepley MPI_Comm comm; 113870034214SMatthew G. Knepley PetscInt nvtxs = numVertices; /* The number of vertices in full graph */ 113970034214SMatthew G. Knepley PetscInt *vtxdist; /* Distribution of vertices across processes */ 114070034214SMatthew G. Knepley PetscInt *xadj = start; /* Start of edge list for each vertex */ 114170034214SMatthew G. Knepley PetscInt *adjncy = adjacency; /* Edge lists for all vertices */ 114270034214SMatthew G. Knepley PetscInt *vwgt = NULL; /* Vertex weights */ 114370034214SMatthew G. Knepley PetscInt *adjwgt = NULL; /* Edge weights */ 114470034214SMatthew G. Knepley PetscInt wgtflag = 0; /* Indicates which weights are present */ 114570034214SMatthew G. Knepley PetscInt numflag = 0; /* Indicates initial offset (0 or 1) */ 114670034214SMatthew G. Knepley PetscInt ncon = 1; /* The number of weights per vertex */ 114770034214SMatthew G. Knepley PetscReal *tpwgts; /* The fraction of vertex weights assigned to each partition */ 114870034214SMatthew G. Knepley PetscReal *ubvec; /* The balance intolerance for vertex weights */ 114970034214SMatthew G. Knepley PetscInt options[5]; /* Options */ 115070034214SMatthew G. Knepley /* Outputs */ 115170034214SMatthew G. Knepley PetscInt edgeCut; /* The number of edges cut by the partition */ 115270034214SMatthew G. Knepley PetscInt *assignment, *points; 115377623264SMatthew G. Knepley PetscMPIInt rank, p, v, i; 115470034214SMatthew G. Knepley PetscErrorCode ierr; 115570034214SMatthew G. Knepley 115670034214SMatthew G. Knepley PetscFunctionBegin; 115777623264SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr); 115870034214SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 115970034214SMatthew G. Knepley options[0] = 0; /* Use all defaults */ 116070034214SMatthew G. Knepley /* Calculate vertex distribution */ 116170034214SMatthew G. Knepley ierr = PetscMalloc4(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment);CHKERRQ(ierr); 116270034214SMatthew G. Knepley vtxdist[0] = 0; 116370034214SMatthew G. Knepley ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr); 116470034214SMatthew G. Knepley for (p = 2; p <= nparts; ++p) { 116570034214SMatthew G. Knepley vtxdist[p] += vtxdist[p-1]; 116670034214SMatthew G. Knepley } 116770034214SMatthew G. Knepley /* Calculate weights */ 116870034214SMatthew G. Knepley for (p = 0; p < nparts; ++p) { 116970034214SMatthew G. Knepley tpwgts[p] = 1.0/nparts; 117070034214SMatthew G. Knepley } 117170034214SMatthew G. Knepley ubvec[0] = 1.05; 117270034214SMatthew G. Knepley 117370034214SMatthew G. Knepley if (nparts == 1) { 117470034214SMatthew G. Knepley ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt)); 117570034214SMatthew G. Knepley } else { 117670034214SMatthew G. Knepley if (vtxdist[1] == vtxdist[nparts]) { 117770034214SMatthew G. Knepley if (!rank) { 117870034214SMatthew G. Knepley PetscStackPush("METIS_PartGraphKway"); 117970034214SMatthew G. Knepley ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment); 118070034214SMatthew G. Knepley PetscStackPop; 118170034214SMatthew G. Knepley if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()"); 118270034214SMatthew G. Knepley } 118370034214SMatthew G. Knepley } else { 118470034214SMatthew G. Knepley PetscStackPush("ParMETIS_V3_PartKway"); 118570034214SMatthew G. Knepley ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm); 118670034214SMatthew G. Knepley PetscStackPop; 118770034214SMatthew G. Knepley if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in ParMETIS_V3_PartKway()"); 118870034214SMatthew G. Knepley } 118970034214SMatthew G. Knepley } 119070034214SMatthew G. Knepley /* Convert to PetscSection+IS */ 119177623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 119277623264SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);} 119377623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 119470034214SMatthew G. Knepley ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr); 119577623264SMatthew G. Knepley for (p = 0, i = 0; p < nparts; ++p) { 119670034214SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) { 119770034214SMatthew G. Knepley if (assignment[v] == p) points[i++] = v; 119870034214SMatthew G. Knepley } 119970034214SMatthew G. Knepley } 120070034214SMatthew G. Knepley if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs); 120170034214SMatthew G. Knepley ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 120270034214SMatthew G. Knepley ierr = PetscFree4(vtxdist,tpwgts,ubvec,assignment);CHKERRQ(ierr); 12039b80ac48SMatthew G. Knepley PetscFunctionReturn(0); 120470034214SMatthew G. Knepley #else 120577623264SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis."); 120670034214SMatthew G. Knepley #endif 120770034214SMatthew G. Knepley } 120870034214SMatthew G. Knepley 120977623264SMatthew G. Knepley #undef __FUNCT__ 121077623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_ParMetis" 121177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part) 121277623264SMatthew G. Knepley { 121377623264SMatthew G. Knepley PetscFunctionBegin; 121477623264SMatthew G. Knepley part->ops->view = PetscPartitionerView_ParMetis; 121577623264SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_ParMetis; 121677623264SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_ParMetis; 121777623264SMatthew G. Knepley PetscFunctionReturn(0); 121877623264SMatthew G. Knepley } 121977623264SMatthew G. Knepley 122077623264SMatthew G. Knepley /*MC 122177623264SMatthew G. Knepley PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMetis library 122277623264SMatthew G. Knepley 122377623264SMatthew G. Knepley Level: intermediate 122477623264SMatthew G. Knepley 122577623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 122677623264SMatthew G. Knepley M*/ 122777623264SMatthew G. Knepley 122877623264SMatthew G. Knepley #undef __FUNCT__ 122977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_ParMetis" 123077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part) 123177623264SMatthew G. Knepley { 123277623264SMatthew G. Knepley PetscPartitioner_ParMetis *p; 123377623264SMatthew G. Knepley PetscErrorCode ierr; 123477623264SMatthew G. Knepley 123577623264SMatthew G. Knepley PetscFunctionBegin; 123677623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 123777623264SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 123877623264SMatthew G. Knepley part->data = p; 123977623264SMatthew G. Knepley 124077623264SMatthew G. Knepley ierr = PetscPartitionerInitialize_ParMetis(part);CHKERRQ(ierr); 124177623264SMatthew G. Knepley ierr = PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionercite);CHKERRQ(ierr); 124270034214SMatthew G. Knepley PetscFunctionReturn(0); 124370034214SMatthew G. Knepley } 124470034214SMatthew G. Knepley 124570034214SMatthew G. Knepley #undef __FUNCT__ 12465680f57bSMatthew G. Knepley #define __FUNCT__ "DMPlexGetPartitioner" 12475680f57bSMatthew G. Knepley /*@ 12485680f57bSMatthew G. Knepley DMPlexGetPartitioner - Get the mesh partitioner 12495680f57bSMatthew G. Knepley 12505680f57bSMatthew G. Knepley Not collective 12515680f57bSMatthew G. Knepley 12525680f57bSMatthew G. Knepley Input Parameter: 12535680f57bSMatthew G. Knepley . dm - The DM 12545680f57bSMatthew G. Knepley 12555680f57bSMatthew G. Knepley Output Parameter: 12565680f57bSMatthew G. Knepley . part - The PetscPartitioner 12575680f57bSMatthew G. Knepley 12585680f57bSMatthew G. Knepley Level: developer 12595680f57bSMatthew G. Knepley 126098599a47SLawrence Mitchell Note: This gets a borrowed reference, so the user should not destroy this PetscPartitioner. 126198599a47SLawrence Mitchell 126298599a47SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexSetPartitioner(), PetscPartitionerCreate() 12635680f57bSMatthew G. Knepley @*/ 12645680f57bSMatthew G. Knepley PetscErrorCode DMPlexGetPartitioner(DM dm, PetscPartitioner *part) 12655680f57bSMatthew G. Knepley { 12665680f57bSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 12675680f57bSMatthew G. Knepley 12685680f57bSMatthew G. Knepley PetscFunctionBegin; 12695680f57bSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12705680f57bSMatthew G. Knepley PetscValidPointer(part, 2); 12715680f57bSMatthew G. Knepley *part = mesh->partitioner; 12725680f57bSMatthew G. Knepley PetscFunctionReturn(0); 12735680f57bSMatthew G. Knepley } 12745680f57bSMatthew G. Knepley 12755680f57bSMatthew G. Knepley #undef __FUNCT__ 127671bb2955SLawrence Mitchell #define __FUNCT__ "DMPlexSetPartitioner" 127771bb2955SLawrence Mitchell /*@ 127871bb2955SLawrence Mitchell DMPlexSetPartitioner - Set the mesh partitioner 127971bb2955SLawrence Mitchell 128071bb2955SLawrence Mitchell logically collective on dm and part 128171bb2955SLawrence Mitchell 128271bb2955SLawrence Mitchell Input Parameters: 128371bb2955SLawrence Mitchell + dm - The DM 128471bb2955SLawrence Mitchell - part - The partitioner 128571bb2955SLawrence Mitchell 128671bb2955SLawrence Mitchell Level: developer 128771bb2955SLawrence Mitchell 128871bb2955SLawrence Mitchell Note: Any existing PetscPartitioner will be destroyed. 128971bb2955SLawrence Mitchell 129071bb2955SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexGetPartitioner(), PetscPartitionerCreate() 129171bb2955SLawrence Mitchell @*/ 129271bb2955SLawrence Mitchell PetscErrorCode DMPlexSetPartitioner(DM dm, PetscPartitioner part) 129371bb2955SLawrence Mitchell { 129471bb2955SLawrence Mitchell DM_Plex *mesh = (DM_Plex *) dm->data; 129571bb2955SLawrence Mitchell PetscErrorCode ierr; 129671bb2955SLawrence Mitchell 129771bb2955SLawrence Mitchell PetscFunctionBegin; 129871bb2955SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 129971bb2955SLawrence Mitchell PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 2); 130071bb2955SLawrence Mitchell ierr = PetscObjectReference((PetscObject)part);CHKERRQ(ierr); 130171bb2955SLawrence Mitchell ierr = PetscPartitionerDestroy(&mesh->partitioner);CHKERRQ(ierr); 130271bb2955SLawrence Mitchell mesh->partitioner = part; 130371bb2955SLawrence Mitchell PetscFunctionReturn(0); 130471bb2955SLawrence Mitchell } 130571bb2955SLawrence Mitchell 130671bb2955SLawrence Mitchell #undef __FUNCT__ 1307270bba0cSToby Isaac #define __FUNCT__ "DMPlexPartitionLabelClosure_Tree" 1308270bba0cSToby Isaac static PetscErrorCode DMPlexPartitionLabelClosure_Tree(DM dm, DMLabel label, PetscInt rank, PetscInt point) 1309270bba0cSToby Isaac { 1310270bba0cSToby Isaac PetscInt parent, closureSize, *closure = NULL, i, pStart, pEnd; 1311270bba0cSToby Isaac PetscErrorCode ierr; 1312270bba0cSToby Isaac 1313270bba0cSToby Isaac PetscFunctionBegin; 1314270bba0cSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr); 1315270bba0cSToby Isaac if (parent == point) PetscFunctionReturn(0); 1316270bba0cSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 1317270bba0cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1318270bba0cSToby Isaac for (i = 0; i < closureSize; i++) { 1319270bba0cSToby Isaac PetscInt cpoint = closure[2*i]; 1320270bba0cSToby Isaac 1321270bba0cSToby Isaac ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr); 1322270bba0cSToby Isaac ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint);CHKERRQ(ierr); 1323270bba0cSToby Isaac } 1324270bba0cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1325270bba0cSToby Isaac PetscFunctionReturn(0); 1326270bba0cSToby Isaac } 1327270bba0cSToby Isaac 1328270bba0cSToby Isaac #undef __FUNCT__ 13295abbe4feSMichael Lange #define __FUNCT__ "DMPlexPartitionLabelClosure" 13305abbe4feSMichael Lange /*@ 13315abbe4feSMichael Lange DMPlexPartitionLabelClosure - Add the closure of all points to the partition label 13325abbe4feSMichael Lange 13335abbe4feSMichael Lange Input Parameters: 13345abbe4feSMichael Lange + dm - The DM 13355abbe4feSMichael Lange - label - DMLabel assinging ranks to remote roots 13365abbe4feSMichael Lange 13375abbe4feSMichael Lange Level: developer 13385abbe4feSMichael Lange 13395abbe4feSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 13405abbe4feSMichael Lange @*/ 13415abbe4feSMichael Lange PetscErrorCode DMPlexPartitionLabelClosure(DM dm, DMLabel label) 13429ffc88e4SToby Isaac { 13435abbe4feSMichael Lange IS rankIS, pointIS; 13445abbe4feSMichael Lange const PetscInt *ranks, *points; 13455abbe4feSMichael Lange PetscInt numRanks, numPoints, r, p, c, closureSize; 13465abbe4feSMichael Lange PetscInt *closure = NULL; 13479ffc88e4SToby Isaac PetscErrorCode ierr; 13489ffc88e4SToby Isaac 13499ffc88e4SToby Isaac PetscFunctionBegin; 13505abbe4feSMichael Lange ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr); 13515abbe4feSMichael Lange ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 13525abbe4feSMichael Lange ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 13535abbe4feSMichael Lange for (r = 0; r < numRanks; ++r) { 13545abbe4feSMichael Lange const PetscInt rank = ranks[r]; 13559ffc88e4SToby Isaac 13565abbe4feSMichael Lange ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr); 13575abbe4feSMichael Lange ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 13585abbe4feSMichael Lange ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 13595abbe4feSMichael Lange for (p = 0; p < numPoints; ++p) { 13605abbe4feSMichael Lange ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1361270bba0cSToby Isaac for (c = 0; c < closureSize*2; c += 2) { 1362270bba0cSToby Isaac ierr = DMLabelSetValue(label, closure[c], rank);CHKERRQ(ierr); 1363270bba0cSToby Isaac ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,closure[c]);CHKERRQ(ierr); 1364270bba0cSToby Isaac } 13659ffc88e4SToby Isaac } 13665abbe4feSMichael Lange ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 13675abbe4feSMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 13689ffc88e4SToby Isaac } 13697de78196SMichael Lange if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);} 13705abbe4feSMichael Lange ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 13715abbe4feSMichael Lange ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 13729ffc88e4SToby Isaac PetscFunctionReturn(0); 13739ffc88e4SToby Isaac } 13749ffc88e4SToby Isaac 13759ffc88e4SToby Isaac #undef __FUNCT__ 137624d039d7SMichael Lange #define __FUNCT__ "DMPlexPartitionLabelAdjacency" 137724d039d7SMichael Lange /*@ 137824d039d7SMichael Lange DMPlexPartitionLabelAdjacency - Add one level of adjacent points to the partition label 137924d039d7SMichael Lange 138024d039d7SMichael Lange Input Parameters: 138124d039d7SMichael Lange + dm - The DM 138224d039d7SMichael Lange - label - DMLabel assinging ranks to remote roots 138324d039d7SMichael Lange 138424d039d7SMichael Lange Level: developer 138524d039d7SMichael Lange 138624d039d7SMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 138724d039d7SMichael Lange @*/ 138824d039d7SMichael Lange PetscErrorCode DMPlexPartitionLabelAdjacency(DM dm, DMLabel label) 138970034214SMatthew G. Knepley { 139024d039d7SMichael Lange IS rankIS, pointIS; 139124d039d7SMichael Lange const PetscInt *ranks, *points; 139224d039d7SMichael Lange PetscInt numRanks, numPoints, r, p, a, adjSize; 139324d039d7SMichael Lange PetscInt *adj = NULL; 139470034214SMatthew G. Knepley PetscErrorCode ierr; 139570034214SMatthew G. Knepley 139670034214SMatthew G. Knepley PetscFunctionBegin; 139724d039d7SMichael Lange ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr); 139824d039d7SMichael Lange ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 139924d039d7SMichael Lange ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 140024d039d7SMichael Lange for (r = 0; r < numRanks; ++r) { 140124d039d7SMichael Lange const PetscInt rank = ranks[r]; 140270034214SMatthew G. Knepley 140324d039d7SMichael Lange ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr); 140424d039d7SMichael Lange ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 140524d039d7SMichael Lange ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 140670034214SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 140724d039d7SMichael Lange adjSize = PETSC_DETERMINE; 140824d039d7SMichael Lange ierr = DMPlexGetAdjacency(dm, points[p], &adjSize, &adj);CHKERRQ(ierr); 140924d039d7SMichael Lange for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(label, adj[a], rank);CHKERRQ(ierr);} 141070034214SMatthew G. Knepley } 141124d039d7SMichael Lange ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 141224d039d7SMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 141370034214SMatthew G. Knepley } 141424d039d7SMichael Lange ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 141524d039d7SMichael Lange ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 141624d039d7SMichael Lange ierr = PetscFree(adj);CHKERRQ(ierr); 141724d039d7SMichael Lange PetscFunctionReturn(0); 141870034214SMatthew G. Knepley } 141970034214SMatthew G. Knepley 142024d039d7SMichael Lange #undef __FUNCT__ 14211fd9873aSMichael Lange #define __FUNCT__ "DMPlexPartitionLabelInvert" 14221fd9873aSMichael Lange /*@ 14231fd9873aSMichael Lange DMPlexPartitionLabelInvert - Create a partition label of remote roots from a local root label 14241fd9873aSMichael Lange 14251fd9873aSMichael Lange Input Parameters: 14261fd9873aSMichael Lange + dm - The DM 14271fd9873aSMichael Lange . rootLabel - DMLabel assinging ranks to local roots 14281fd9873aSMichael Lange . processSF - A star forest mapping into the local index on each remote rank 14291fd9873aSMichael Lange 14301fd9873aSMichael Lange Output Parameter: 14311fd9873aSMichael Lange - leafLabel - DMLabel assinging ranks to remote roots 14321fd9873aSMichael Lange 14331fd9873aSMichael Lange Note: The rootLabel defines a send pattern by mapping local points to remote target ranks. The 14341fd9873aSMichael Lange resulting leafLabel is a receiver mapping of remote roots to their parent rank. 14351fd9873aSMichael Lange 14361fd9873aSMichael Lange Level: developer 14371fd9873aSMichael Lange 14381fd9873aSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 14391fd9873aSMichael Lange @*/ 14401fd9873aSMichael Lange PetscErrorCode DMPlexPartitionLabelInvert(DM dm, DMLabel rootLabel, PetscSF processSF, DMLabel leafLabel) 14411fd9873aSMichael Lange { 14421fd9873aSMichael Lange MPI_Comm comm; 14431fd9873aSMichael Lange PetscMPIInt rank, numProcs; 14441fd9873aSMichael Lange PetscInt p, n, numNeighbors, size, l, nleaves; 14451fd9873aSMichael Lange PetscSF sfPoint; 14461fd9873aSMichael Lange PetscSFNode *rootPoints, *leafPoints; 14471fd9873aSMichael Lange PetscSection rootSection, leafSection; 14481fd9873aSMichael Lange const PetscSFNode *remote; 14491fd9873aSMichael Lange const PetscInt *local, *neighbors; 14501fd9873aSMichael Lange IS valueIS; 14511fd9873aSMichael Lange PetscErrorCode ierr; 14521fd9873aSMichael Lange 14531fd9873aSMichael Lange PetscFunctionBegin; 14541fd9873aSMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 14551fd9873aSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 14561fd9873aSMichael Lange ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 14571fd9873aSMichael Lange ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 14581fd9873aSMichael Lange 14591fd9873aSMichael Lange /* Convert to (point, rank) and use actual owners */ 14601fd9873aSMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 14611fd9873aSMichael Lange ierr = PetscSectionSetChart(rootSection, 0, numProcs);CHKERRQ(ierr); 14621fd9873aSMichael Lange ierr = DMLabelGetValueIS(rootLabel, &valueIS);CHKERRQ(ierr); 14631fd9873aSMichael Lange ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr); 14641fd9873aSMichael Lange ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr); 14651fd9873aSMichael Lange for (n = 0; n < numNeighbors; ++n) { 14661fd9873aSMichael Lange PetscInt numPoints; 14671fd9873aSMichael Lange 14681fd9873aSMichael Lange ierr = DMLabelGetStratumSize(rootLabel, neighbors[n], &numPoints);CHKERRQ(ierr); 14691fd9873aSMichael Lange ierr = PetscSectionAddDof(rootSection, neighbors[n], numPoints);CHKERRQ(ierr); 14701fd9873aSMichael Lange } 14711fd9873aSMichael Lange ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 14721fd9873aSMichael Lange ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr); 14731fd9873aSMichael Lange ierr = PetscMalloc1(size, &rootPoints);CHKERRQ(ierr); 14741fd9873aSMichael Lange ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 14751fd9873aSMichael Lange for (n = 0; n < numNeighbors; ++n) { 14761fd9873aSMichael Lange IS pointIS; 14771fd9873aSMichael Lange const PetscInt *points; 14781fd9873aSMichael Lange PetscInt off, numPoints, p; 14791fd9873aSMichael Lange 14801fd9873aSMichael Lange ierr = PetscSectionGetOffset(rootSection, neighbors[n], &off);CHKERRQ(ierr); 14811fd9873aSMichael Lange ierr = DMLabelGetStratumIS(rootLabel, neighbors[n], &pointIS);CHKERRQ(ierr); 14821fd9873aSMichael Lange ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 14831fd9873aSMichael Lange ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 14841fd9873aSMichael Lange for (p = 0; p < numPoints; ++p) { 1485f8987ae8SMichael Lange if (local) {ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);} 1486f8987ae8SMichael Lange else {l = -1;} 14871fd9873aSMichael Lange if (l >= 0) {rootPoints[off+p] = remote[l];} 14881fd9873aSMichael Lange else {rootPoints[off+p].index = points[p]; rootPoints[off+p].rank = rank;} 14891fd9873aSMichael Lange } 14901fd9873aSMichael Lange ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 14911fd9873aSMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 14921fd9873aSMichael Lange } 14931fd9873aSMichael Lange ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr); 14941fd9873aSMichael Lange ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 14951fd9873aSMichael Lange /* Communicate overlap */ 14961fd9873aSMichael Lange ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); 14971fd9873aSMichael Lange ierr = DMPlexDistributeData(dm, processSF, rootSection, MPIU_2INT, rootPoints, leafSection, (void**) &leafPoints);CHKERRQ(ierr); 14981fd9873aSMichael Lange /* Filter remote contributions (ovLeafPoints) into the overlapSF */ 14991fd9873aSMichael Lange ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr); 15001fd9873aSMichael Lange for (p = 0; p < size; p++) { 15011fd9873aSMichael Lange ierr = DMLabelSetValue(leafLabel, leafPoints[p].index, leafPoints[p].rank);CHKERRQ(ierr); 15021fd9873aSMichael Lange } 15031fd9873aSMichael Lange ierr = PetscFree(rootPoints);CHKERRQ(ierr); 15041fd9873aSMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 15051fd9873aSMichael Lange ierr = PetscFree(leafPoints);CHKERRQ(ierr); 15061fd9873aSMichael Lange ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 15071fd9873aSMichael Lange PetscFunctionReturn(0); 15081fd9873aSMichael Lange } 15091fd9873aSMichael Lange 15101fd9873aSMichael Lange #undef __FUNCT__ 1511aa3148a8SMichael Lange #define __FUNCT__ "DMPlexPartitionLabelCreateSF" 1512aa3148a8SMichael Lange /*@ 1513aa3148a8SMichael Lange DMPlexPartitionLabelCreateSF - Create a star forest from a label that assigns ranks to points 1514aa3148a8SMichael Lange 1515aa3148a8SMichael Lange Input Parameters: 1516aa3148a8SMichael Lange + dm - The DM 1517aa3148a8SMichael Lange . label - DMLabel assinging ranks to remote roots 1518aa3148a8SMichael Lange 1519aa3148a8SMichael Lange Output Parameter: 1520aa3148a8SMichael Lange - sf - The star forest communication context encapsulating the defined mapping 1521aa3148a8SMichael Lange 1522aa3148a8SMichael Lange Note: The incoming label is a receiver mapping of remote points to their parent rank. 1523aa3148a8SMichael Lange 1524aa3148a8SMichael Lange Level: developer 1525aa3148a8SMichael Lange 1526aa3148a8SMichael Lange .seealso: DMPlexDistribute(), DMPlexCreateOverlap 1527aa3148a8SMichael Lange @*/ 1528aa3148a8SMichael Lange PetscErrorCode DMPlexPartitionLabelCreateSF(DM dm, DMLabel label, PetscSF *sf) 1529aa3148a8SMichael Lange { 153043f7d02bSMichael Lange PetscMPIInt rank, numProcs; 153143f7d02bSMichael Lange PetscInt n, numRemote, p, numPoints, pStart, pEnd, idx = 0; 1532aa3148a8SMichael Lange PetscSFNode *remotePoints; 153343f7d02bSMichael Lange IS remoteRootIS; 153443f7d02bSMichael Lange const PetscInt *remoteRoots; 1535aa3148a8SMichael Lange PetscErrorCode ierr; 1536aa3148a8SMichael Lange 1537aa3148a8SMichael Lange PetscFunctionBegin; 153843f7d02bSMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1539aa3148a8SMichael Lange ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 1540aa3148a8SMichael Lange 1541aa3148a8SMichael Lange for (numRemote = 0, n = 0; n < numProcs; ++n) { 1542aa3148a8SMichael Lange ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr); 1543aa3148a8SMichael Lange numRemote += numPoints; 1544aa3148a8SMichael Lange } 1545aa3148a8SMichael Lange ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr); 154643f7d02bSMichael Lange /* Put owned points first */ 154743f7d02bSMichael Lange ierr = DMLabelGetStratumSize(label, rank, &numPoints);CHKERRQ(ierr); 154843f7d02bSMichael Lange if (numPoints > 0) { 154943f7d02bSMichael Lange ierr = DMLabelGetStratumIS(label, rank, &remoteRootIS);CHKERRQ(ierr); 155043f7d02bSMichael Lange ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 155143f7d02bSMichael Lange for (p = 0; p < numPoints; p++) { 155243f7d02bSMichael Lange remotePoints[idx].index = remoteRoots[p]; 155343f7d02bSMichael Lange remotePoints[idx].rank = rank; 155443f7d02bSMichael Lange idx++; 155543f7d02bSMichael Lange } 155643f7d02bSMichael Lange ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 155743f7d02bSMichael Lange ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr); 155843f7d02bSMichael Lange } 155943f7d02bSMichael Lange /* Now add remote points */ 156043f7d02bSMichael Lange for (n = 0; n < numProcs; ++n) { 1561aa3148a8SMichael Lange ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr); 156243f7d02bSMichael Lange if (numPoints <= 0 || n == rank) continue; 1563aa3148a8SMichael Lange ierr = DMLabelGetStratumIS(label, n, &remoteRootIS);CHKERRQ(ierr); 1564aa3148a8SMichael Lange ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 1565aa3148a8SMichael Lange for (p = 0; p < numPoints; p++) { 1566aa3148a8SMichael Lange remotePoints[idx].index = remoteRoots[p]; 1567aa3148a8SMichael Lange remotePoints[idx].rank = n; 1568aa3148a8SMichael Lange idx++; 1569aa3148a8SMichael Lange } 1570aa3148a8SMichael Lange ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 1571aa3148a8SMichael Lange ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr); 1572aa3148a8SMichael Lange } 1573aa3148a8SMichael Lange ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), sf);CHKERRQ(ierr); 1574aa3148a8SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1575aa3148a8SMichael Lange ierr = PetscSFSetGraph(*sf, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 157670034214SMatthew G. Knepley PetscFunctionReturn(0); 157770034214SMatthew G. Knepley } 1578