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 29532c4e7dSMichael Lange /*@C 30532c4e7dSMichael Lange DMPlexCreatePartitionerGraph - Create a CSR graph of point connections for the partitioner 31532c4e7dSMichael Lange 32532c4e7dSMichael Lange Input Parameters: 33532c4e7dSMichael Lange + dm - The mesh DM dm 34532c4e7dSMichael Lange - height - Height of the strata from which to construct the graph 35532c4e7dSMichael Lange 36532c4e7dSMichael Lange Output Parameter: 37532c4e7dSMichael Lange + numVertices - Number of vertices in the graph 383fa7752dSToby Isaac . offsets - Point offsets in the graph 393fa7752dSToby Isaac . adjacency - Point connectivity in the graph 403fa7752dSToby Isaac - globalNumbering - A map from the local cell numbering to the global numbering used in "adjacency". Negative indicates that the cell is a duplicate from another process. 41532c4e7dSMichael Lange 42532c4e7dSMichael Lange The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 43532c4e7dSMichael Lange DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 44532c4e7dSMichael Lange representation on the mesh. 45532c4e7dSMichael Lange 46532c4e7dSMichael Lange Level: developer 47532c4e7dSMichael Lange 48532c4e7dSMichael Lange .seealso: PetscPartitionerGetType(), PetscPartitionerCreate(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 49532c4e7dSMichael Lange @*/ 503fa7752dSToby Isaac PetscErrorCode DMPlexCreatePartitionerGraph(DM dm, PetscInt height, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency, IS *globalNumbering) 51532c4e7dSMichael Lange { 5243f7d02bSMichael Lange PetscInt p, pStart, pEnd, a, adjSize, idx, size, nroots; 53389e55d8SMichael Lange PetscInt *adj = NULL, *vOffsets = NULL, *graph = NULL; 548cfe4c1fSMichael Lange IS cellNumbering; 558cfe4c1fSMichael Lange const PetscInt *cellNum; 56532c4e7dSMichael Lange PetscBool useCone, useClosure; 57532c4e7dSMichael Lange PetscSection section; 58532c4e7dSMichael Lange PetscSegBuffer adjBuffer; 598cfe4c1fSMichael Lange PetscSF sfPoint; 60532c4e7dSMichael Lange PetscMPIInt rank; 61532c4e7dSMichael Lange PetscErrorCode ierr; 62532c4e7dSMichael Lange 63532c4e7dSMichael Lange PetscFunctionBegin; 64532c4e7dSMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 65532c4e7dSMichael Lange ierr = DMPlexGetHeightStratum(dm, height, &pStart, &pEnd);CHKERRQ(ierr); 668cfe4c1fSMichael Lange ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 678cfe4c1fSMichael Lange ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 68532c4e7dSMichael Lange /* Build adjacency graph via a section/segbuffer */ 69532c4e7dSMichael Lange ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 70532c4e7dSMichael Lange ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 71532c4e7dSMichael Lange ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&adjBuffer);CHKERRQ(ierr); 72532c4e7dSMichael Lange /* Always use FVM adjacency to create partitioner graph */ 73532c4e7dSMichael Lange ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr); 74532c4e7dSMichael Lange ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr); 75532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr); 76532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);CHKERRQ(ierr); 77f0927f4eSMatthew G. Knepley ierr = DMPlexCreateCellNumbering_Internal(dm, PETSC_TRUE, &cellNumbering);CHKERRQ(ierr); 783fa7752dSToby Isaac if (globalNumbering) { 793fa7752dSToby Isaac ierr = PetscObjectReference((PetscObject)cellNumbering);CHKERRQ(ierr); 803fa7752dSToby Isaac *globalNumbering = cellNumbering; 813fa7752dSToby Isaac } 828cfe4c1fSMichael Lange ierr = ISGetIndices(cellNumbering, &cellNum);CHKERRQ(ierr); 838cfe4c1fSMichael Lange for (*numVertices = 0, p = pStart; p < pEnd; p++) { 848cfe4c1fSMichael Lange /* Skip non-owned cells in parallel (ParMetis expects no overlap) */ 858cfe4c1fSMichael Lange if (nroots > 0) {if (cellNum[p] < 0) continue;} 86532c4e7dSMichael Lange adjSize = PETSC_DETERMINE; 87532c4e7dSMichael Lange ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 88532c4e7dSMichael Lange for (a = 0; a < adjSize; ++a) { 89532c4e7dSMichael Lange const PetscInt point = adj[a]; 90532c4e7dSMichael Lange if (point != p && pStart <= point && point < pEnd) { 91532c4e7dSMichael Lange PetscInt *PETSC_RESTRICT pBuf; 92532c4e7dSMichael Lange ierr = PetscSectionAddDof(section, p, 1);CHKERRQ(ierr); 93532c4e7dSMichael Lange ierr = PetscSegBufferGetInts(adjBuffer, 1, &pBuf);CHKERRQ(ierr); 94532c4e7dSMichael Lange *pBuf = point; 95532c4e7dSMichael Lange } 96532c4e7dSMichael Lange } 978cfe4c1fSMichael Lange (*numVertices)++; 98532c4e7dSMichael Lange } 99532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr); 100532c4e7dSMichael Lange ierr = DMPlexSetAdjacencyUseClosure(dm, useClosure);CHKERRQ(ierr); 101532c4e7dSMichael Lange /* Derive CSR graph from section/segbuffer */ 102532c4e7dSMichael Lange ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 103532c4e7dSMichael Lange ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr); 104389e55d8SMichael Lange ierr = PetscMalloc1(*numVertices+1, &vOffsets);CHKERRQ(ierr); 10543f7d02bSMichael Lange for (idx = 0, p = pStart; p < pEnd; p++) { 10643f7d02bSMichael Lange if (nroots > 0) {if (cellNum[p] < 0) continue;} 10743f7d02bSMichael Lange ierr = PetscSectionGetOffset(section, p, &(vOffsets[idx++]));CHKERRQ(ierr); 10843f7d02bSMichael Lange } 109532c4e7dSMichael Lange vOffsets[*numVertices] = size; 110532c4e7dSMichael Lange if (offsets) *offsets = vOffsets; 111389e55d8SMichael Lange ierr = PetscSegBufferExtractAlloc(adjBuffer, &graph);CHKERRQ(ierr); 112bf4602e4SToby Isaac { 1138cfe4c1fSMichael Lange ISLocalToGlobalMapping ltogCells; 1148cfe4c1fSMichael Lange PetscInt n, size, *cells_arr; 1158cfe4c1fSMichael Lange /* In parallel, apply a global cell numbering to the graph */ 1168cfe4c1fSMichael Lange ierr = ISRestoreIndices(cellNumbering, &cellNum);CHKERRQ(ierr); 1178cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingCreateIS(cellNumbering, <ogCells);CHKERRQ(ierr); 1188cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingGetSize(ltogCells, &size);CHKERRQ(ierr); 1198cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingGetIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr); 1208cfe4c1fSMichael Lange /* Convert to positive global cell numbers */ 1218cfe4c1fSMichael Lange for (n=0; n<size; n++) {if (cells_arr[n] < 0) cells_arr[n] = -(cells_arr[n]+1);} 1228cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingRestoreIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr); 123389e55d8SMichael Lange ierr = ISLocalToGlobalMappingApplyBlock(ltogCells, vOffsets[*numVertices], graph, graph);CHKERRQ(ierr); 1248cfe4c1fSMichael Lange ierr = ISLocalToGlobalMappingDestroy(<ogCells);CHKERRQ(ierr); 125f0927f4eSMatthew G. Knepley ierr = ISDestroy(&cellNumbering);CHKERRQ(ierr); 1268cfe4c1fSMichael Lange } 127389e55d8SMichael Lange if (adjacency) *adjacency = graph; 128532c4e7dSMichael Lange /* Clean up */ 129532c4e7dSMichael Lange ierr = PetscSegBufferDestroy(&adjBuffer);CHKERRQ(ierr); 130532c4e7dSMichael Lange ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 131532c4e7dSMichael Lange ierr = PetscFree(adj);CHKERRQ(ierr); 132532c4e7dSMichael Lange PetscFunctionReturn(0); 133532c4e7dSMichael Lange } 134532c4e7dSMichael Lange 135d5577e40SMatthew G. Knepley /*@C 136d5577e40SMatthew G. Knepley DMPlexCreateNeighborCSR - Create a mesh graph (cell-cell adjacency) in parallel CSR format. 137d5577e40SMatthew G. Knepley 138d5577e40SMatthew G. Knepley Collective 139d5577e40SMatthew G. Knepley 140d5577e40SMatthew G. Knepley Input Arguments: 141d5577e40SMatthew G. Knepley + dm - The DMPlex 142d5577e40SMatthew G. Knepley - cellHeight - The height of mesh points to treat as cells (default should be 0) 143d5577e40SMatthew G. Knepley 144d5577e40SMatthew G. Knepley Output Arguments: 145d5577e40SMatthew G. Knepley + numVertices - The number of local vertices in the graph, or cells in the mesh. 146d5577e40SMatthew G. Knepley . offsets - The offset to the adjacency list for each cell 147d5577e40SMatthew G. Knepley - adjacency - The adjacency list for all cells 148d5577e40SMatthew G. Knepley 149d5577e40SMatthew G. Knepley Note: This is suitable for input to a mesh partitioner like ParMetis. 150d5577e40SMatthew G. Knepley 151d5577e40SMatthew G. Knepley Level: advanced 152d5577e40SMatthew G. Knepley 153d5577e40SMatthew G. Knepley .seealso: DMPlexCreate() 154d5577e40SMatthew G. Knepley @*/ 15570034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency) 15670034214SMatthew G. Knepley { 15770034214SMatthew G. Knepley const PetscInt maxFaceCases = 30; 15870034214SMatthew G. Knepley PetscInt numFaceCases = 0; 15970034214SMatthew G. Knepley PetscInt numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */ 16070034214SMatthew G. Knepley PetscInt *off, *adj; 16170034214SMatthew G. Knepley PetscInt *neighborCells = NULL; 16270034214SMatthew G. Knepley PetscInt dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell; 16370034214SMatthew G. Knepley PetscErrorCode ierr; 16470034214SMatthew G. Knepley 16570034214SMatthew G. Knepley PetscFunctionBegin; 16670034214SMatthew G. Knepley /* For parallel partitioning, I think you have to communicate supports */ 167c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 16870034214SMatthew G. Knepley cellDim = dim - cellHeight; 16970034214SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 17070034214SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr); 17170034214SMatthew G. Knepley if (cEnd - cStart == 0) { 17270034214SMatthew G. Knepley if (numVertices) *numVertices = 0; 17370034214SMatthew G. Knepley if (offsets) *offsets = NULL; 17470034214SMatthew G. Knepley if (adjacency) *adjacency = NULL; 17570034214SMatthew G. Knepley PetscFunctionReturn(0); 17670034214SMatthew G. Knepley } 17770034214SMatthew G. Knepley numCells = cEnd - cStart; 17870034214SMatthew G. Knepley faceDepth = depth - cellHeight; 17970034214SMatthew G. Knepley if (dim == depth) { 18070034214SMatthew G. Knepley PetscInt f, fStart, fEnd; 18170034214SMatthew G. Knepley 18270034214SMatthew G. Knepley ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr); 18370034214SMatthew G. Knepley /* Count neighboring cells */ 18470034214SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr); 18570034214SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 18670034214SMatthew G. Knepley const PetscInt *support; 18770034214SMatthew G. Knepley PetscInt supportSize; 18870034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 18970034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr); 19070034214SMatthew G. Knepley if (supportSize == 2) { 1919ffc88e4SToby Isaac PetscInt numChildren; 1929ffc88e4SToby Isaac 1939ffc88e4SToby Isaac ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 1949ffc88e4SToby Isaac if (!numChildren) { 19570034214SMatthew G. Knepley ++off[support[0]-cStart+1]; 19670034214SMatthew G. Knepley ++off[support[1]-cStart+1]; 19770034214SMatthew G. Knepley } 19870034214SMatthew G. Knepley } 1999ffc88e4SToby Isaac } 20070034214SMatthew G. Knepley /* Prefix sum */ 20170034214SMatthew G. Knepley for (c = 1; c <= numCells; ++c) off[c] += off[c-1]; 20270034214SMatthew G. Knepley if (adjacency) { 20370034214SMatthew G. Knepley PetscInt *tmp; 20470034214SMatthew G. Knepley 20570034214SMatthew G. Knepley ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr); 206854ce69bSBarry Smith ierr = PetscMalloc1(numCells+1, &tmp);CHKERRQ(ierr); 20770034214SMatthew G. Knepley ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr); 20870034214SMatthew G. Knepley /* Get neighboring cells */ 20970034214SMatthew G. Knepley for (f = fStart; f < fEnd; ++f) { 21070034214SMatthew G. Knepley const PetscInt *support; 21170034214SMatthew G. Knepley PetscInt supportSize; 21270034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr); 21370034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr); 21470034214SMatthew G. Knepley if (supportSize == 2) { 2159ffc88e4SToby Isaac PetscInt numChildren; 2169ffc88e4SToby Isaac 2179ffc88e4SToby Isaac ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr); 2189ffc88e4SToby Isaac if (!numChildren) { 21970034214SMatthew G. Knepley adj[tmp[support[0]-cStart]++] = support[1]; 22070034214SMatthew G. Knepley adj[tmp[support[1]-cStart]++] = support[0]; 22170034214SMatthew G. Knepley } 22270034214SMatthew G. Knepley } 2239ffc88e4SToby Isaac } 22470034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG) 22570034214SMatthew 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); 22670034214SMatthew G. Knepley #endif 22770034214SMatthew G. Knepley ierr = PetscFree(tmp);CHKERRQ(ierr); 22870034214SMatthew G. Knepley } 22970034214SMatthew G. Knepley if (numVertices) *numVertices = numCells; 23070034214SMatthew G. Knepley if (offsets) *offsets = off; 23170034214SMatthew G. Knepley if (adjacency) *adjacency = adj; 23270034214SMatthew G. Knepley PetscFunctionReturn(0); 23370034214SMatthew G. Knepley } 23470034214SMatthew G. Knepley /* Setup face recognition */ 23570034214SMatthew G. Knepley if (faceDepth == 1) { 23670034214SMatthew 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 */ 23770034214SMatthew G. Knepley 23870034214SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 23970034214SMatthew G. Knepley PetscInt corners; 24070034214SMatthew G. Knepley 24170034214SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr); 24270034214SMatthew G. Knepley if (!cornersSeen[corners]) { 24370034214SMatthew G. Knepley PetscInt nFV; 24470034214SMatthew G. Knepley 24570034214SMatthew G. Knepley if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases"); 24670034214SMatthew G. Knepley cornersSeen[corners] = 1; 24770034214SMatthew G. Knepley 24870034214SMatthew G. Knepley ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr); 24970034214SMatthew G. Knepley 25070034214SMatthew G. Knepley numFaceVertices[numFaceCases++] = nFV; 25170034214SMatthew G. Knepley } 25270034214SMatthew G. Knepley } 25370034214SMatthew G. Knepley } 25470034214SMatthew G. Knepley ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr); 25570034214SMatthew G. Knepley /* Count neighboring cells */ 25670034214SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 25770034214SMatthew G. Knepley PetscInt numNeighbors = PETSC_DETERMINE, n; 25870034214SMatthew G. Knepley 2598b0b4c70SToby Isaac ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr); 26070034214SMatthew G. Knepley /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */ 26170034214SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 26270034214SMatthew G. Knepley PetscInt cellPair[2]; 26370034214SMatthew G. Knepley PetscBool found = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE; 26470034214SMatthew G. Knepley PetscInt meetSize = 0; 26570034214SMatthew G. Knepley const PetscInt *meet = NULL; 26670034214SMatthew G. Knepley 26770034214SMatthew G. Knepley cellPair[0] = cell; cellPair[1] = neighborCells[n]; 26870034214SMatthew G. Knepley if (cellPair[0] == cellPair[1]) continue; 26970034214SMatthew G. Knepley if (!found) { 27070034214SMatthew G. Knepley ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 27170034214SMatthew G. Knepley if (meetSize) { 27270034214SMatthew G. Knepley PetscInt f; 27370034214SMatthew G. Knepley 27470034214SMatthew G. Knepley for (f = 0; f < numFaceCases; ++f) { 27570034214SMatthew G. Knepley if (numFaceVertices[f] == meetSize) { 27670034214SMatthew G. Knepley found = PETSC_TRUE; 27770034214SMatthew G. Knepley break; 27870034214SMatthew G. Knepley } 27970034214SMatthew G. Knepley } 28070034214SMatthew G. Knepley } 28170034214SMatthew G. Knepley ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 28270034214SMatthew G. Knepley } 28370034214SMatthew G. Knepley if (found) ++off[cell-cStart+1]; 28470034214SMatthew G. Knepley } 28570034214SMatthew G. Knepley } 28670034214SMatthew G. Knepley /* Prefix sum */ 28770034214SMatthew G. Knepley for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1]; 28870034214SMatthew G. Knepley 28970034214SMatthew G. Knepley if (adjacency) { 29070034214SMatthew G. Knepley ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr); 29170034214SMatthew G. Knepley /* Get neighboring cells */ 29270034214SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 29370034214SMatthew G. Knepley PetscInt numNeighbors = PETSC_DETERMINE, n; 29470034214SMatthew G. Knepley PetscInt cellOffset = 0; 29570034214SMatthew G. Knepley 2968b0b4c70SToby Isaac ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr); 29770034214SMatthew G. Knepley /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */ 29870034214SMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 29970034214SMatthew G. Knepley PetscInt cellPair[2]; 30070034214SMatthew G. Knepley PetscBool found = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE; 30170034214SMatthew G. Knepley PetscInt meetSize = 0; 30270034214SMatthew G. Knepley const PetscInt *meet = NULL; 30370034214SMatthew G. Knepley 30470034214SMatthew G. Knepley cellPair[0] = cell; cellPair[1] = neighborCells[n]; 30570034214SMatthew G. Knepley if (cellPair[0] == cellPair[1]) continue; 30670034214SMatthew G. Knepley if (!found) { 30770034214SMatthew G. Knepley ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 30870034214SMatthew G. Knepley if (meetSize) { 30970034214SMatthew G. Knepley PetscInt f; 31070034214SMatthew G. Knepley 31170034214SMatthew G. Knepley for (f = 0; f < numFaceCases; ++f) { 31270034214SMatthew G. Knepley if (numFaceVertices[f] == meetSize) { 31370034214SMatthew G. Knepley found = PETSC_TRUE; 31470034214SMatthew G. Knepley break; 31570034214SMatthew G. Knepley } 31670034214SMatthew G. Knepley } 31770034214SMatthew G. Knepley } 31870034214SMatthew G. Knepley ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr); 31970034214SMatthew G. Knepley } 32070034214SMatthew G. Knepley if (found) { 32170034214SMatthew G. Knepley adj[off[cell-cStart]+cellOffset] = neighborCells[n]; 32270034214SMatthew G. Knepley ++cellOffset; 32370034214SMatthew G. Knepley } 32470034214SMatthew G. Knepley } 32570034214SMatthew G. Knepley } 32670034214SMatthew G. Knepley } 32770034214SMatthew G. Knepley ierr = PetscFree(neighborCells);CHKERRQ(ierr); 32870034214SMatthew G. Knepley if (numVertices) *numVertices = numCells; 32970034214SMatthew G. Knepley if (offsets) *offsets = off; 33070034214SMatthew G. Knepley if (adjacency) *adjacency = adj; 33170034214SMatthew G. Knepley PetscFunctionReturn(0); 33270034214SMatthew G. Knepley } 33370034214SMatthew G. Knepley 33477623264SMatthew G. Knepley /*@C 33577623264SMatthew G. Knepley PetscPartitionerRegister - Adds a new PetscPartitioner implementation 33677623264SMatthew G. Knepley 33777623264SMatthew G. Knepley Not Collective 33877623264SMatthew G. Knepley 33977623264SMatthew G. Knepley Input Parameters: 34077623264SMatthew G. Knepley + name - The name of a new user-defined creation routine 34177623264SMatthew G. Knepley - create_func - The creation routine itself 34277623264SMatthew G. Knepley 34377623264SMatthew G. Knepley Notes: 34477623264SMatthew G. Knepley PetscPartitionerRegister() may be called multiple times to add several user-defined PetscPartitioners 34577623264SMatthew G. Knepley 34677623264SMatthew G. Knepley Sample usage: 34777623264SMatthew G. Knepley .vb 34877623264SMatthew G. Knepley PetscPartitionerRegister("my_part", MyPetscPartitionerCreate); 34977623264SMatthew G. Knepley .ve 35077623264SMatthew G. Knepley 35177623264SMatthew G. Knepley Then, your PetscPartitioner type can be chosen with the procedural interface via 35277623264SMatthew G. Knepley .vb 35377623264SMatthew G. Knepley PetscPartitionerCreate(MPI_Comm, PetscPartitioner *); 35477623264SMatthew G. Knepley PetscPartitionerSetType(PetscPartitioner, "my_part"); 35577623264SMatthew G. Knepley .ve 35677623264SMatthew G. Knepley or at runtime via the option 35777623264SMatthew G. Knepley .vb 35877623264SMatthew G. Knepley -petscpartitioner_type my_part 35977623264SMatthew G. Knepley .ve 36077623264SMatthew G. Knepley 36177623264SMatthew G. Knepley Level: advanced 36277623264SMatthew G. Knepley 36377623264SMatthew G. Knepley .keywords: PetscPartitioner, register 36477623264SMatthew G. Knepley .seealso: PetscPartitionerRegisterAll(), PetscPartitionerRegisterDestroy() 36577623264SMatthew G. Knepley 36677623264SMatthew G. Knepley @*/ 36777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerRegister(const char sname[], PetscErrorCode (*function)(PetscPartitioner)) 36877623264SMatthew G. Knepley { 36977623264SMatthew G. Knepley PetscErrorCode ierr; 37077623264SMatthew G. Knepley 37177623264SMatthew G. Knepley PetscFunctionBegin; 37277623264SMatthew G. Knepley ierr = PetscFunctionListAdd(&PetscPartitionerList, sname, function);CHKERRQ(ierr); 37377623264SMatthew G. Knepley PetscFunctionReturn(0); 37477623264SMatthew G. Knepley } 37577623264SMatthew G. Knepley 37677623264SMatthew G. Knepley /*@C 37777623264SMatthew G. Knepley PetscPartitionerSetType - Builds a particular PetscPartitioner 37877623264SMatthew G. Knepley 37977623264SMatthew G. Knepley Collective on PetscPartitioner 38077623264SMatthew G. Knepley 38177623264SMatthew G. Knepley Input Parameters: 38277623264SMatthew G. Knepley + part - The PetscPartitioner object 38377623264SMatthew G. Knepley - name - The kind of partitioner 38477623264SMatthew G. Knepley 38577623264SMatthew G. Knepley Options Database Key: 38677623264SMatthew G. Knepley . -petscpartitioner_type <type> - Sets the PetscPartitioner type; use -help for a list of available types 38777623264SMatthew G. Knepley 38877623264SMatthew G. Knepley Level: intermediate 38977623264SMatthew G. Knepley 39077623264SMatthew G. Knepley .keywords: PetscPartitioner, set, type 39177623264SMatthew G. Knepley .seealso: PetscPartitionerGetType(), PetscPartitionerCreate() 39277623264SMatthew G. Knepley @*/ 39377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name) 39477623264SMatthew G. Knepley { 39577623264SMatthew G. Knepley PetscErrorCode (*r)(PetscPartitioner); 39677623264SMatthew G. Knepley PetscBool match; 39777623264SMatthew G. Knepley PetscErrorCode ierr; 39877623264SMatthew G. Knepley 39977623264SMatthew G. Knepley PetscFunctionBegin; 40077623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 40177623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) part, name, &match);CHKERRQ(ierr); 40277623264SMatthew G. Knepley if (match) PetscFunctionReturn(0); 40377623264SMatthew G. Knepley 4040f51fdf8SToby Isaac ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr); 40577623264SMatthew G. Knepley ierr = PetscFunctionListFind(PetscPartitionerList, name, &r);CHKERRQ(ierr); 40677623264SMatthew G. Knepley if (!r) SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name); 40777623264SMatthew G. Knepley 40877623264SMatthew G. Knepley if (part->ops->destroy) { 40977623264SMatthew G. Knepley ierr = (*part->ops->destroy)(part);CHKERRQ(ierr); 41077623264SMatthew G. Knepley part->ops->destroy = NULL; 41177623264SMatthew G. Knepley } 41277623264SMatthew G. Knepley ierr = (*r)(part);CHKERRQ(ierr); 41377623264SMatthew G. Knepley ierr = PetscObjectChangeTypeName((PetscObject) part, name);CHKERRQ(ierr); 41477623264SMatthew G. Knepley PetscFunctionReturn(0); 41577623264SMatthew G. Knepley } 41677623264SMatthew G. Knepley 41777623264SMatthew G. Knepley /*@C 41877623264SMatthew G. Knepley PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object. 41977623264SMatthew G. Knepley 42077623264SMatthew G. Knepley Not Collective 42177623264SMatthew G. Knepley 42277623264SMatthew G. Knepley Input Parameter: 42377623264SMatthew G. Knepley . part - The PetscPartitioner 42477623264SMatthew G. Knepley 42577623264SMatthew G. Knepley Output Parameter: 42677623264SMatthew G. Knepley . name - The PetscPartitioner type name 42777623264SMatthew G. Knepley 42877623264SMatthew G. Knepley Level: intermediate 42977623264SMatthew G. Knepley 43077623264SMatthew G. Knepley .keywords: PetscPartitioner, get, type, name 43177623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitionerCreate() 43277623264SMatthew G. Knepley @*/ 43377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name) 43477623264SMatthew G. Knepley { 43577623264SMatthew G. Knepley PetscErrorCode ierr; 43677623264SMatthew G. Knepley 43777623264SMatthew G. Knepley PetscFunctionBegin; 43877623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 43977623264SMatthew G. Knepley PetscValidPointer(name, 2); 4400f51fdf8SToby Isaac ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr); 44177623264SMatthew G. Knepley *name = ((PetscObject) part)->type_name; 44277623264SMatthew G. Knepley PetscFunctionReturn(0); 44377623264SMatthew G. Knepley } 44477623264SMatthew G. Knepley 44577623264SMatthew G. Knepley /*@C 44677623264SMatthew G. Knepley PetscPartitionerView - Views a PetscPartitioner 44777623264SMatthew G. Knepley 44877623264SMatthew G. Knepley Collective on PetscPartitioner 44977623264SMatthew G. Knepley 45077623264SMatthew G. Knepley Input Parameter: 45177623264SMatthew G. Knepley + part - the PetscPartitioner object to view 45277623264SMatthew G. Knepley - v - the viewer 45377623264SMatthew G. Knepley 45477623264SMatthew G. Knepley Level: developer 45577623264SMatthew G. Knepley 45677623264SMatthew G. Knepley .seealso: PetscPartitionerDestroy() 45777623264SMatthew G. Knepley @*/ 45877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v) 45977623264SMatthew G. Knepley { 46077623264SMatthew G. Knepley PetscErrorCode ierr; 46177623264SMatthew G. Knepley 46277623264SMatthew G. Knepley PetscFunctionBegin; 46377623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 46477623264SMatthew G. Knepley if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) part), &v);CHKERRQ(ierr);} 46577623264SMatthew G. Knepley if (part->ops->view) {ierr = (*part->ops->view)(part, v);CHKERRQ(ierr);} 46677623264SMatthew G. Knepley PetscFunctionReturn(0); 46777623264SMatthew G. Knepley } 46877623264SMatthew G. Knepley 46977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner part) 47077623264SMatthew G. Knepley { 47177623264SMatthew G. Knepley const char *defaultType; 47277623264SMatthew G. Knepley char name[256]; 47377623264SMatthew G. Knepley PetscBool flg; 47477623264SMatthew G. Knepley PetscErrorCode ierr; 47577623264SMatthew G. Knepley 47677623264SMatthew G. Knepley PetscFunctionBegin; 47777623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 4787cf71cfdSLisandro Dalcin if (!((PetscObject) part)->type_name) 4797cf71cfdSLisandro Dalcin #if defined(PETSC_HAVE_CHACO) 4807cf71cfdSLisandro Dalcin defaultType = PETSCPARTITIONERCHACO; 4817cf71cfdSLisandro Dalcin #elif defined(PETSC_HAVE_PARMETIS) 4827cf71cfdSLisandro Dalcin defaultType = PETSCPARTITIONERPARMETIS; 4837cf71cfdSLisandro Dalcin #else 4847cf71cfdSLisandro Dalcin defaultType = PETSCPARTITIONERSIMPLE; 4857cf71cfdSLisandro Dalcin #endif 4867cf71cfdSLisandro Dalcin else 4877cf71cfdSLisandro Dalcin defaultType = ((PetscObject) part)->type_name; 4880f51fdf8SToby Isaac ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr); 48977623264SMatthew G. Knepley 49077623264SMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr); 49177623264SMatthew G. Knepley ierr = PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, defaultType, name, 256, &flg);CHKERRQ(ierr); 49277623264SMatthew G. Knepley if (flg) { 49377623264SMatthew G. Knepley ierr = PetscPartitionerSetType(part, name);CHKERRQ(ierr); 49477623264SMatthew G. Knepley } else if (!((PetscObject) part)->type_name) { 49577623264SMatthew G. Knepley ierr = PetscPartitionerSetType(part, defaultType);CHKERRQ(ierr); 49677623264SMatthew G. Knepley } 49777623264SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 49877623264SMatthew G. Knepley PetscFunctionReturn(0); 49977623264SMatthew G. Knepley } 50077623264SMatthew G. Knepley 50177623264SMatthew G. Knepley /*@ 50277623264SMatthew G. Knepley PetscPartitionerSetFromOptions - sets parameters in a PetscPartitioner from the options database 50377623264SMatthew G. Knepley 50477623264SMatthew G. Knepley Collective on PetscPartitioner 50577623264SMatthew G. Knepley 50677623264SMatthew G. Knepley Input Parameter: 50777623264SMatthew G. Knepley . part - the PetscPartitioner object to set options for 50877623264SMatthew G. Knepley 50977623264SMatthew G. Knepley Level: developer 51077623264SMatthew G. Knepley 51177623264SMatthew G. Knepley .seealso: PetscPartitionerView() 51277623264SMatthew G. Knepley @*/ 51377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part) 51477623264SMatthew G. Knepley { 51577623264SMatthew G. Knepley PetscErrorCode ierr; 51677623264SMatthew G. Knepley 51777623264SMatthew G. Knepley PetscFunctionBegin; 51877623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 51977623264SMatthew G. Knepley ierr = PetscPartitionerSetTypeFromOptions_Internal(part);CHKERRQ(ierr); 52077623264SMatthew G. Knepley 52177623264SMatthew G. Knepley ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr); 522077101c0SMatthew G. Knepley if (part->ops->setfromoptions) {ierr = (*part->ops->setfromoptions)(PetscOptionsObject,part);CHKERRQ(ierr);} 52377623264SMatthew G. Knepley /* process any options handlers added with PetscObjectAddOptionsHandler() */ 5240633abcbSJed Brown ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) part);CHKERRQ(ierr); 52577623264SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 52677623264SMatthew G. Knepley ierr = PetscPartitionerViewFromOptions(part, NULL, "-petscpartitioner_view");CHKERRQ(ierr); 52777623264SMatthew G. Knepley PetscFunctionReturn(0); 52877623264SMatthew G. Knepley } 52977623264SMatthew G. Knepley 53077623264SMatthew G. Knepley /*@C 53177623264SMatthew G. Knepley PetscPartitionerSetUp - Construct data structures for the PetscPartitioner 53277623264SMatthew G. Knepley 53377623264SMatthew G. Knepley Collective on PetscPartitioner 53477623264SMatthew G. Knepley 53577623264SMatthew G. Knepley Input Parameter: 53677623264SMatthew G. Knepley . part - the PetscPartitioner object to setup 53777623264SMatthew G. Knepley 53877623264SMatthew G. Knepley Level: developer 53977623264SMatthew G. Knepley 54077623264SMatthew G. Knepley .seealso: PetscPartitionerView(), PetscPartitionerDestroy() 54177623264SMatthew G. Knepley @*/ 54277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part) 54377623264SMatthew G. Knepley { 54477623264SMatthew G. Knepley PetscErrorCode ierr; 54577623264SMatthew G. Knepley 54677623264SMatthew G. Knepley PetscFunctionBegin; 54777623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 54877623264SMatthew G. Knepley if (part->ops->setup) {ierr = (*part->ops->setup)(part);CHKERRQ(ierr);} 54977623264SMatthew G. Knepley PetscFunctionReturn(0); 55077623264SMatthew G. Knepley } 55177623264SMatthew G. Knepley 55277623264SMatthew G. Knepley /*@ 55377623264SMatthew G. Knepley PetscPartitionerDestroy - Destroys a PetscPartitioner object 55477623264SMatthew G. Knepley 55577623264SMatthew G. Knepley Collective on PetscPartitioner 55677623264SMatthew G. Knepley 55777623264SMatthew G. Knepley Input Parameter: 55877623264SMatthew G. Knepley . part - the PetscPartitioner object to destroy 55977623264SMatthew G. Knepley 56077623264SMatthew G. Knepley Level: developer 56177623264SMatthew G. Knepley 56277623264SMatthew G. Knepley .seealso: PetscPartitionerView() 56377623264SMatthew G. Knepley @*/ 56477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part) 56577623264SMatthew G. Knepley { 56677623264SMatthew G. Knepley PetscErrorCode ierr; 56777623264SMatthew G. Knepley 56877623264SMatthew G. Knepley PetscFunctionBegin; 56977623264SMatthew G. Knepley if (!*part) PetscFunctionReturn(0); 57077623264SMatthew G. Knepley PetscValidHeaderSpecific((*part), PETSCPARTITIONER_CLASSID, 1); 57177623264SMatthew G. Knepley 57277623264SMatthew G. Knepley if (--((PetscObject)(*part))->refct > 0) {*part = 0; PetscFunctionReturn(0);} 57377623264SMatthew G. Knepley ((PetscObject) (*part))->refct = 0; 57477623264SMatthew G. Knepley 57577623264SMatthew G. Knepley if ((*part)->ops->destroy) {ierr = (*(*part)->ops->destroy)(*part);CHKERRQ(ierr);} 57677623264SMatthew G. Knepley ierr = PetscHeaderDestroy(part);CHKERRQ(ierr); 57777623264SMatthew G. Knepley PetscFunctionReturn(0); 57877623264SMatthew G. Knepley } 57977623264SMatthew G. Knepley 58077623264SMatthew G. Knepley /*@ 58177623264SMatthew G. Knepley PetscPartitionerCreate - Creates an empty PetscPartitioner object. The type can then be set with PetscPartitionerSetType(). 58277623264SMatthew G. Knepley 58377623264SMatthew G. Knepley Collective on MPI_Comm 58477623264SMatthew G. Knepley 58577623264SMatthew G. Knepley Input Parameter: 58677623264SMatthew G. Knepley . comm - The communicator for the PetscPartitioner object 58777623264SMatthew G. Knepley 58877623264SMatthew G. Knepley Output Parameter: 58977623264SMatthew G. Knepley . part - The PetscPartitioner object 59077623264SMatthew G. Knepley 59177623264SMatthew G. Knepley Level: beginner 59277623264SMatthew G. Knepley 593dae52e14SToby Isaac .seealso: PetscPartitionerSetType(), PETSCPARTITIONERCHACO, PETSCPARTITIONERPARMETIS, PETSCPARTITIONERSHELL, PETSCPARTITIONERSIMPLE, PETSCPARTITIONERGATHER 59477623264SMatthew G. Knepley @*/ 59577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part) 59677623264SMatthew G. Knepley { 59777623264SMatthew G. Knepley PetscPartitioner p; 59877623264SMatthew G. Knepley PetscErrorCode ierr; 59977623264SMatthew G. Knepley 60077623264SMatthew G. Knepley PetscFunctionBegin; 60177623264SMatthew G. Knepley PetscValidPointer(part, 2); 60277623264SMatthew G. Knepley *part = NULL; 60383cde681SMatthew G. Knepley ierr = DMInitializePackage();CHKERRQ(ierr); 60477623264SMatthew G. Knepley 60573107ff1SLisandro Dalcin ierr = PetscHeaderCreate(p, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView);CHKERRQ(ierr); 60677623264SMatthew G. Knepley 60777623264SMatthew G. Knepley *part = p; 60877623264SMatthew G. Knepley PetscFunctionReturn(0); 60977623264SMatthew G. Knepley } 61077623264SMatthew G. Knepley 61177623264SMatthew G. Knepley /*@ 61277623264SMatthew G. Knepley PetscPartitionerPartition - Create a non-overlapping partition of the cells in the mesh 61377623264SMatthew G. Knepley 61477623264SMatthew G. Knepley Collective on DM 61577623264SMatthew G. Knepley 61677623264SMatthew G. Knepley Input Parameters: 61777623264SMatthew G. Knepley + part - The PetscPartitioner 618f8987ae8SMichael Lange - dm - The mesh DM 61977623264SMatthew G. Knepley 62077623264SMatthew G. Knepley Output Parameters: 62177623264SMatthew G. Knepley + partSection - The PetscSection giving the division of points by partition 622f8987ae8SMichael Lange - partition - The list of points by partition 62377623264SMatthew G. Knepley 62477623264SMatthew G. Knepley Note: Instead of cells, points at a given height can be partitioned by calling PetscPartitionerSetPointHeight() 62577623264SMatthew G. Knepley 62677623264SMatthew G. Knepley Level: developer 62777623264SMatthew G. Knepley 62877623264SMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerSetPointHeight(), PetscPartitionerCreate() 6294b15ede2SMatthew G. Knepley @*/ 630f8987ae8SMichael Lange PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, DM dm, PetscSection partSection, IS *partition) 63177623264SMatthew G. Knepley { 63277623264SMatthew G. Knepley PetscMPIInt size; 63377623264SMatthew G. Knepley PetscErrorCode ierr; 63477623264SMatthew G. Knepley 63577623264SMatthew G. Knepley PetscFunctionBegin; 63677623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 63777623264SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 63877623264SMatthew G. Knepley PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 4); 63977623264SMatthew G. Knepley PetscValidPointer(partition, 5); 64077623264SMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) part), &size);CHKERRQ(ierr); 64177623264SMatthew G. Knepley if (size == 1) { 64277623264SMatthew G. Knepley PetscInt *points; 64377623264SMatthew G. Knepley PetscInt cStart, cEnd, c; 64477623264SMatthew G. Knepley 64577623264SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr); 64677623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, size);CHKERRQ(ierr); 64777623264SMatthew G. Knepley ierr = PetscSectionSetDof(partSection, 0, cEnd-cStart);CHKERRQ(ierr); 64877623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 64977623264SMatthew G. Knepley ierr = PetscMalloc1(cEnd-cStart, &points);CHKERRQ(ierr); 65077623264SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) points[c] = c; 65177623264SMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 65277623264SMatthew G. Knepley PetscFunctionReturn(0); 65377623264SMatthew G. Knepley } 654*13fef8f3SMatthew G. Knepley /* Obvious cheating, but I cannot think of a better way to do this. The DMSetFromOptions() has refinement in it, 655*13fef8f3SMatthew G. Knepley so we cannot call it and have it feed down to the partitioner before partitioning */ 656*13fef8f3SMatthew G. Knepley ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 65777623264SMatthew G. Knepley if (part->height == 0) { 65877623264SMatthew G. Knepley PetscInt numVertices; 65977623264SMatthew G. Knepley PetscInt *start = NULL; 66077623264SMatthew G. Knepley PetscInt *adjacency = NULL; 6613fa7752dSToby Isaac IS globalNumbering; 66277623264SMatthew G. Knepley 6633fa7752dSToby Isaac ierr = DMPlexCreatePartitionerGraph(dm, 0, &numVertices, &start, &adjacency, &globalNumbering);CHKERRQ(ierr); 66477623264SMatthew G. Knepley if (!part->ops->partition) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONGSTATE, "PetscPartitioner has no type"); 66577623264SMatthew G. Knepley ierr = (*part->ops->partition)(part, dm, size, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr); 66677623264SMatthew G. Knepley ierr = PetscFree(start);CHKERRQ(ierr); 66777623264SMatthew G. Knepley ierr = PetscFree(adjacency);CHKERRQ(ierr); 6683fa7752dSToby Isaac if (globalNumbering) { /* partition is wrt global unique numbering: change this to be wrt local numbering */ 6693fa7752dSToby Isaac const PetscInt *globalNum; 6703fa7752dSToby Isaac const PetscInt *partIdx; 6713fa7752dSToby Isaac PetscInt *map, cStart, cEnd; 6723fa7752dSToby Isaac PetscInt *adjusted, i, localSize, offset; 6733fa7752dSToby Isaac IS newPartition; 6743fa7752dSToby Isaac 6753fa7752dSToby Isaac ierr = ISGetLocalSize(*partition,&localSize);CHKERRQ(ierr); 6763fa7752dSToby Isaac ierr = PetscMalloc1(localSize,&adjusted);CHKERRQ(ierr); 6773fa7752dSToby Isaac ierr = ISGetIndices(globalNumbering,&globalNum);CHKERRQ(ierr); 6783fa7752dSToby Isaac ierr = ISGetIndices(*partition,&partIdx);CHKERRQ(ierr); 6793fa7752dSToby Isaac ierr = PetscMalloc1(localSize,&map);CHKERRQ(ierr); 6803fa7752dSToby Isaac ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr); 6813fa7752dSToby Isaac for (i = cStart, offset = 0; i < cEnd; i++) { 6823fa7752dSToby Isaac if (globalNum[i - cStart] >= 0) map[offset++] = i; 6833fa7752dSToby Isaac } 6843fa7752dSToby Isaac for (i = 0; i < localSize; i++) { 6853fa7752dSToby Isaac adjusted[i] = map[partIdx[i]]; 6863fa7752dSToby Isaac } 6873fa7752dSToby Isaac ierr = PetscFree(map);CHKERRQ(ierr); 6883fa7752dSToby Isaac ierr = ISRestoreIndices(*partition,&partIdx);CHKERRQ(ierr); 6893fa7752dSToby Isaac ierr = ISRestoreIndices(globalNumbering,&globalNum);CHKERRQ(ierr); 6903fa7752dSToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,localSize,adjusted,PETSC_OWN_POINTER,&newPartition);CHKERRQ(ierr); 6913fa7752dSToby Isaac ierr = ISDestroy(&globalNumbering);CHKERRQ(ierr); 6923fa7752dSToby Isaac ierr = ISDestroy(partition);CHKERRQ(ierr); 6933fa7752dSToby Isaac *partition = newPartition; 6943fa7752dSToby Isaac } 69577623264SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_OUTOFRANGE, "Invalid height %D for points to partition", part->height); 69677623264SMatthew G. Knepley PetscFunctionReturn(0); 6973fa7752dSToby Isaac 69877623264SMatthew G. Knepley } 69977623264SMatthew G. Knepley 700d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part) 70177623264SMatthew G. Knepley { 70277623264SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 70377623264SMatthew G. Knepley PetscErrorCode ierr; 70477623264SMatthew G. Knepley 70577623264SMatthew G. Knepley PetscFunctionBegin; 70677623264SMatthew G. Knepley ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr); 70777623264SMatthew G. Knepley ierr = ISDestroy(&p->partition);CHKERRQ(ierr); 70877623264SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 70977623264SMatthew G. Knepley PetscFunctionReturn(0); 71077623264SMatthew G. Knepley } 71177623264SMatthew G. Knepley 712d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Shell_Ascii(PetscPartitioner part, PetscViewer viewer) 71377623264SMatthew G. Knepley { 714077101c0SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 71577623264SMatthew G. Knepley PetscViewerFormat format; 71677623264SMatthew G. Knepley PetscErrorCode ierr; 71777623264SMatthew G. Knepley 71877623264SMatthew G. Knepley PetscFunctionBegin; 71977623264SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 72077623264SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Shell Graph Partitioner:\n");CHKERRQ(ierr); 721077101c0SMatthew G. Knepley if (p->random) { 722077101c0SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 723077101c0SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "using random partition\n");CHKERRQ(ierr); 724077101c0SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 725077101c0SMatthew G. Knepley } 72677623264SMatthew G. Knepley PetscFunctionReturn(0); 72777623264SMatthew G. Knepley } 72877623264SMatthew G. Knepley 729d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer) 73077623264SMatthew G. Knepley { 73177623264SMatthew G. Knepley PetscBool iascii; 73277623264SMatthew G. Knepley PetscErrorCode ierr; 73377623264SMatthew G. Knepley 73477623264SMatthew G. Knepley PetscFunctionBegin; 73577623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 73677623264SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 73777623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 73877623264SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_Shell_Ascii(part, viewer);CHKERRQ(ierr);} 73977623264SMatthew G. Knepley PetscFunctionReturn(0); 74077623264SMatthew G. Knepley } 74177623264SMatthew G. Knepley 742d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscOptionItems *PetscOptionsObject, PetscPartitioner part) 743077101c0SMatthew G. Knepley { 744077101c0SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 745077101c0SMatthew G. Knepley PetscErrorCode ierr; 746077101c0SMatthew G. Knepley 747077101c0SMatthew G. Knepley PetscFunctionBegin; 748077101c0SMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject, "PetscPartitionerShell Options");CHKERRQ(ierr); 749077101c0SMatthew G. Knepley ierr = PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &p->random, NULL);CHKERRQ(ierr); 750077101c0SMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 751077101c0SMatthew G. Knepley PetscFunctionReturn(0); 752077101c0SMatthew G. Knepley } 753077101c0SMatthew G. Knepley 754d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 75577623264SMatthew G. Knepley { 75677623264SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 75777623264SMatthew G. Knepley PetscInt np; 75877623264SMatthew G. Knepley PetscErrorCode ierr; 75977623264SMatthew G. Knepley 76077623264SMatthew G. Knepley PetscFunctionBegin; 761077101c0SMatthew G. Knepley if (p->random) { 762077101c0SMatthew G. Knepley PetscRandom r; 763077101c0SMatthew G. Knepley PetscInt *sizes, *points, v; 764077101c0SMatthew G. Knepley 765077101c0SMatthew G. Knepley ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 766c717d290SMatthew G. Knepley ierr = PetscRandomSetInterval(r, 0.0, (PetscScalar) nparts);CHKERRQ(ierr); 767077101c0SMatthew G. Knepley ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 768077101c0SMatthew G. Knepley ierr = PetscCalloc2(nparts, &sizes, numVertices, &points);CHKERRQ(ierr); 769077101c0SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 770077101c0SMatthew G. Knepley PetscReal val; 771077101c0SMatthew G. Knepley PetscInt part; 772077101c0SMatthew G. Knepley 773077101c0SMatthew G. Knepley ierr = PetscRandomGetValueReal(r, &val);CHKERRQ(ierr); 774077101c0SMatthew G. Knepley part = PetscFloorReal(val); 775077101c0SMatthew G. Knepley points[v] = part; 776077101c0SMatthew G. Knepley ++sizes[part]; 777077101c0SMatthew G. Knepley } 778077101c0SMatthew G. Knepley ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 779077101c0SMatthew G. Knepley ierr = PetscPartitionerShellSetPartition(part, nparts, sizes, points);CHKERRQ(ierr); 780077101c0SMatthew G. Knepley ierr = PetscFree2(sizes, points);CHKERRQ(ierr); 781077101c0SMatthew G. Knepley } 782c717d290SMatthew G. Knepley if (!p->section) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()"); 78377623264SMatthew G. Knepley ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr); 78477623264SMatthew G. Knepley if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np); 78577623264SMatthew G. Knepley ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr); 78677623264SMatthew G. Knepley if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np); 7875680f57bSMatthew G. Knepley ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr); 78877623264SMatthew G. Knepley *partition = p->partition; 78977623264SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr); 79077623264SMatthew G. Knepley PetscFunctionReturn(0); 79177623264SMatthew G. Knepley } 79277623264SMatthew G. Knepley 793d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part) 79477623264SMatthew G. Knepley { 79577623264SMatthew G. Knepley PetscFunctionBegin; 79677623264SMatthew G. Knepley part->ops->view = PetscPartitionerView_Shell; 797077101c0SMatthew G. Knepley part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell; 79877623264SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_Shell; 79977623264SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_Shell; 80077623264SMatthew G. Knepley PetscFunctionReturn(0); 80177623264SMatthew G. Knepley } 80277623264SMatthew G. Knepley 80377623264SMatthew G. Knepley /*MC 80477623264SMatthew G. Knepley PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object 80577623264SMatthew G. Knepley 80677623264SMatthew G. Knepley Level: intermediate 80777623264SMatthew G. Knepley 80877623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 80977623264SMatthew G. Knepley M*/ 81077623264SMatthew G. Knepley 81177623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part) 81277623264SMatthew G. Knepley { 81377623264SMatthew G. Knepley PetscPartitioner_Shell *p; 81477623264SMatthew G. Knepley PetscErrorCode ierr; 81577623264SMatthew G. Knepley 81677623264SMatthew G. Knepley PetscFunctionBegin; 81777623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 81877623264SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 81977623264SMatthew G. Knepley part->data = p; 82077623264SMatthew G. Knepley 82177623264SMatthew G. Knepley ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr); 822077101c0SMatthew G. Knepley p->random = PETSC_FALSE; 82377623264SMatthew G. Knepley PetscFunctionReturn(0); 82477623264SMatthew G. Knepley } 82577623264SMatthew G. Knepley 8265680f57bSMatthew G. Knepley /*@C 8275680f57bSMatthew G. Knepley PetscPartitionerShellSetPartition - Set an artifical partition for a mesh 8285680f57bSMatthew G. Knepley 829077101c0SMatthew G. Knepley Collective on Part 8305680f57bSMatthew G. Knepley 8315680f57bSMatthew G. Knepley Input Parameters: 8325680f57bSMatthew G. Knepley + part - The PetscPartitioner 8339852e123SBarry Smith . size - The number of partitions 8349852e123SBarry Smith . sizes - array of size size (or NULL) providing the number of points in each partition 8359758bf69SToby Isaac - points - array of size sum(sizes) (may be NULL iff sizes is NULL), a permutation of the points that groups those assigned to each partition in order (i.e., partition 0 first, partition 1 next, etc.) 8365680f57bSMatthew G. Knepley 8375680f57bSMatthew G. Knepley Level: developer 8385680f57bSMatthew G. Knepley 839b7e49471SLawrence Mitchell Notes: 840b7e49471SLawrence Mitchell 841b7e49471SLawrence Mitchell It is safe to free the sizes and points arrays after use in this routine. 842b7e49471SLawrence Mitchell 8435680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate() 8445680f57bSMatthew G. Knepley @*/ 8459852e123SBarry Smith PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[]) 8465680f57bSMatthew G. Knepley { 8475680f57bSMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 8485680f57bSMatthew G. Knepley PetscInt proc, numPoints; 8495680f57bSMatthew G. Knepley PetscErrorCode ierr; 8505680f57bSMatthew G. Knepley 8515680f57bSMatthew G. Knepley PetscFunctionBegin; 8525680f57bSMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 8535680f57bSMatthew G. Knepley if (sizes) {PetscValidPointer(sizes, 3);} 854c717d290SMatthew G. Knepley if (points) {PetscValidPointer(points, 4);} 8555680f57bSMatthew G. Knepley ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr); 8565680f57bSMatthew G. Knepley ierr = ISDestroy(&p->partition);CHKERRQ(ierr); 8575680f57bSMatthew G. Knepley ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr); 8589852e123SBarry Smith ierr = PetscSectionSetChart(p->section, 0, size);CHKERRQ(ierr); 8595680f57bSMatthew G. Knepley if (sizes) { 8609852e123SBarry Smith for (proc = 0; proc < size; ++proc) { 8615680f57bSMatthew G. Knepley ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr); 8625680f57bSMatthew G. Knepley } 8635680f57bSMatthew G. Knepley } 8645680f57bSMatthew G. Knepley ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr); 8655680f57bSMatthew G. Knepley ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr); 8665680f57bSMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr); 8675680f57bSMatthew G. Knepley PetscFunctionReturn(0); 8685680f57bSMatthew G. Knepley } 8695680f57bSMatthew G. Knepley 870077101c0SMatthew G. Knepley /*@ 871077101c0SMatthew G. Knepley PetscPartitionerShellSetRandom - Set the flag to use a random partition 872077101c0SMatthew G. Knepley 873077101c0SMatthew G. Knepley Collective on Part 874077101c0SMatthew G. Knepley 875077101c0SMatthew G. Knepley Input Parameters: 876077101c0SMatthew G. Knepley + part - The PetscPartitioner 877077101c0SMatthew G. Knepley - random - The flag to use a random partition 878077101c0SMatthew G. Knepley 879077101c0SMatthew G. Knepley Level: intermediate 880077101c0SMatthew G. Knepley 881077101c0SMatthew G. Knepley .seealso PetscPartitionerShellGetRandom(), PetscPartitionerCreate() 882077101c0SMatthew G. Knepley @*/ 883077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random) 884077101c0SMatthew G. Knepley { 885077101c0SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 886077101c0SMatthew G. Knepley 887077101c0SMatthew G. Knepley PetscFunctionBegin; 888077101c0SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 889077101c0SMatthew G. Knepley p->random = random; 890077101c0SMatthew G. Knepley PetscFunctionReturn(0); 891077101c0SMatthew G. Knepley } 892077101c0SMatthew G. Knepley 893077101c0SMatthew G. Knepley /*@ 894077101c0SMatthew G. Knepley PetscPartitionerShellGetRandom - get the flag to use a random partition 895077101c0SMatthew G. Knepley 896077101c0SMatthew G. Knepley Collective on Part 897077101c0SMatthew G. Knepley 898077101c0SMatthew G. Knepley Input Parameter: 899077101c0SMatthew G. Knepley . part - The PetscPartitioner 900077101c0SMatthew G. Knepley 901077101c0SMatthew G. Knepley Output Parameter 902077101c0SMatthew G. Knepley . random - The flag to use a random partition 903077101c0SMatthew G. Knepley 904077101c0SMatthew G. Knepley Level: intermediate 905077101c0SMatthew G. Knepley 906077101c0SMatthew G. Knepley .seealso PetscPartitionerShellSetRandom(), PetscPartitionerCreate() 907077101c0SMatthew G. Knepley @*/ 908077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random) 909077101c0SMatthew G. Knepley { 910077101c0SMatthew G. Knepley PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data; 911077101c0SMatthew G. Knepley 912077101c0SMatthew G. Knepley PetscFunctionBegin; 913077101c0SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 914077101c0SMatthew G. Knepley PetscValidPointer(random, 2); 915077101c0SMatthew G. Knepley *random = p->random; 916077101c0SMatthew G. Knepley PetscFunctionReturn(0); 917077101c0SMatthew G. Knepley } 918077101c0SMatthew G. Knepley 919d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part) 920555a9cf8SMatthew G. Knepley { 921555a9cf8SMatthew G. Knepley PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 922555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 923555a9cf8SMatthew G. Knepley 924555a9cf8SMatthew G. Knepley PetscFunctionBegin; 925555a9cf8SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 926555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 927555a9cf8SMatthew G. Knepley } 928555a9cf8SMatthew G. Knepley 929d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Simple_Ascii(PetscPartitioner part, PetscViewer viewer) 930555a9cf8SMatthew G. Knepley { 931555a9cf8SMatthew G. Knepley PetscViewerFormat format; 932555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 933555a9cf8SMatthew G. Knepley 934555a9cf8SMatthew G. Knepley PetscFunctionBegin; 935555a9cf8SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 936555a9cf8SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Simple Graph Partitioner:\n");CHKERRQ(ierr); 937555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 938555a9cf8SMatthew G. Knepley } 939555a9cf8SMatthew G. Knepley 940d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer) 941555a9cf8SMatthew G. Knepley { 942555a9cf8SMatthew G. Knepley PetscBool iascii; 943555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 944555a9cf8SMatthew G. Knepley 945555a9cf8SMatthew G. Knepley PetscFunctionBegin; 946555a9cf8SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 947555a9cf8SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 948555a9cf8SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 949555a9cf8SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_Simple_Ascii(part, viewer);CHKERRQ(ierr);} 950555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 951555a9cf8SMatthew G. Knepley } 952555a9cf8SMatthew G. Knepley 953d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 954555a9cf8SMatthew G. Knepley { 955cead94edSToby Isaac MPI_Comm comm; 956555a9cf8SMatthew G. Knepley PetscInt np; 957cead94edSToby Isaac PetscMPIInt size; 958555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 959555a9cf8SMatthew G. Knepley 960555a9cf8SMatthew G. Knepley PetscFunctionBegin; 961cead94edSToby Isaac comm = PetscObjectComm((PetscObject)dm); 962cead94edSToby Isaac ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 963555a9cf8SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 964555a9cf8SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr); 965cead94edSToby Isaac if (size == 1) { 966cead94edSToby Isaac for (np = 0; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np));CHKERRQ(ierr);} 967cead94edSToby Isaac } 968cead94edSToby Isaac else { 969cead94edSToby Isaac PetscMPIInt rank; 970cead94edSToby Isaac PetscInt nvGlobal, *offsets, myFirst, myLast; 971cead94edSToby Isaac 972a679a563SToby Isaac ierr = PetscMalloc1(size+1,&offsets);CHKERRQ(ierr); 973cead94edSToby Isaac offsets[0] = 0; 974cead94edSToby Isaac ierr = MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm);CHKERRQ(ierr); 975cead94edSToby Isaac for (np = 2; np <= size; np++) { 976cead94edSToby Isaac offsets[np] += offsets[np-1]; 977cead94edSToby Isaac } 978cead94edSToby Isaac nvGlobal = offsets[size]; 979cead94edSToby Isaac ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 980cead94edSToby Isaac myFirst = offsets[rank]; 981cead94edSToby Isaac myLast = offsets[rank + 1] - 1; 982cead94edSToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 983cead94edSToby Isaac if (numVertices) { 984cead94edSToby Isaac PetscInt firstPart = 0, firstLargePart = 0; 985cead94edSToby Isaac PetscInt lastPart = 0, lastLargePart = 0; 986cead94edSToby Isaac PetscInt rem = nvGlobal % nparts; 987cead94edSToby Isaac PetscInt pSmall = nvGlobal/nparts; 988cead94edSToby Isaac PetscInt pBig = nvGlobal/nparts + 1; 989cead94edSToby Isaac 990cead94edSToby Isaac 991cead94edSToby Isaac if (rem) { 992cead94edSToby Isaac firstLargePart = myFirst / pBig; 993cead94edSToby Isaac lastLargePart = myLast / pBig; 994cead94edSToby Isaac 995cead94edSToby Isaac if (firstLargePart < rem) { 996cead94edSToby Isaac firstPart = firstLargePart; 997cead94edSToby Isaac } 998cead94edSToby Isaac else { 999cead94edSToby Isaac firstPart = rem + (myFirst - (rem * pBig)) / pSmall; 1000cead94edSToby Isaac } 1001cead94edSToby Isaac if (lastLargePart < rem) { 1002cead94edSToby Isaac lastPart = lastLargePart; 1003cead94edSToby Isaac } 1004cead94edSToby Isaac else { 1005cead94edSToby Isaac lastPart = rem + (myLast - (rem * pBig)) / pSmall; 1006cead94edSToby Isaac } 1007cead94edSToby Isaac } 1008cead94edSToby Isaac else { 1009cead94edSToby Isaac firstPart = myFirst / (nvGlobal/nparts); 1010cead94edSToby Isaac lastPart = myLast / (nvGlobal/nparts); 1011cead94edSToby Isaac } 1012cead94edSToby Isaac 1013cead94edSToby Isaac for (np = firstPart; np <= lastPart; np++) { 1014cead94edSToby Isaac PetscInt PartStart = np * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np); 1015cead94edSToby Isaac PetscInt PartEnd = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1); 1016cead94edSToby Isaac 1017cead94edSToby Isaac PartStart = PetscMax(PartStart,myFirst); 1018cead94edSToby Isaac PartEnd = PetscMin(PartEnd,myLast+1); 1019cead94edSToby Isaac ierr = PetscSectionSetDof(partSection,np,PartEnd-PartStart);CHKERRQ(ierr); 1020cead94edSToby Isaac } 1021cead94edSToby Isaac } 1022cead94edSToby Isaac } 1023cead94edSToby Isaac ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 1024555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 1025555a9cf8SMatthew G. Knepley } 1026555a9cf8SMatthew G. Knepley 1027d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part) 1028555a9cf8SMatthew G. Knepley { 1029555a9cf8SMatthew G. Knepley PetscFunctionBegin; 1030555a9cf8SMatthew G. Knepley part->ops->view = PetscPartitionerView_Simple; 1031555a9cf8SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_Simple; 1032555a9cf8SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_Simple; 1033555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 1034555a9cf8SMatthew G. Knepley } 1035555a9cf8SMatthew G. Knepley 1036555a9cf8SMatthew G. Knepley /*MC 1037555a9cf8SMatthew G. Knepley PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object 1038555a9cf8SMatthew G. Knepley 1039555a9cf8SMatthew G. Knepley Level: intermediate 1040555a9cf8SMatthew G. Knepley 1041555a9cf8SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 1042555a9cf8SMatthew G. Knepley M*/ 1043555a9cf8SMatthew G. Knepley 1044555a9cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part) 1045555a9cf8SMatthew G. Knepley { 1046555a9cf8SMatthew G. Knepley PetscPartitioner_Simple *p; 1047555a9cf8SMatthew G. Knepley PetscErrorCode ierr; 1048555a9cf8SMatthew G. Knepley 1049555a9cf8SMatthew G. Knepley PetscFunctionBegin; 1050555a9cf8SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 1051555a9cf8SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 1052555a9cf8SMatthew G. Knepley part->data = p; 1053555a9cf8SMatthew G. Knepley 1054555a9cf8SMatthew G. Knepley ierr = PetscPartitionerInitialize_Simple(part);CHKERRQ(ierr); 1055555a9cf8SMatthew G. Knepley PetscFunctionReturn(0); 1056555a9cf8SMatthew G. Knepley } 1057555a9cf8SMatthew G. Knepley 1058d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Gather(PetscPartitioner part) 1059dae52e14SToby Isaac { 1060dae52e14SToby Isaac PetscPartitioner_Gather *p = (PetscPartitioner_Gather *) part->data; 1061dae52e14SToby Isaac PetscErrorCode ierr; 1062dae52e14SToby Isaac 1063dae52e14SToby Isaac PetscFunctionBegin; 1064dae52e14SToby Isaac ierr = PetscFree(p);CHKERRQ(ierr); 1065dae52e14SToby Isaac PetscFunctionReturn(0); 1066dae52e14SToby Isaac } 1067dae52e14SToby Isaac 1068d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Gather_Ascii(PetscPartitioner part, PetscViewer viewer) 1069dae52e14SToby Isaac { 1070dae52e14SToby Isaac PetscViewerFormat format; 1071dae52e14SToby Isaac PetscErrorCode ierr; 1072dae52e14SToby Isaac 1073dae52e14SToby Isaac PetscFunctionBegin; 1074dae52e14SToby Isaac ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 1075dae52e14SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Gather Graph Partitioner:\n");CHKERRQ(ierr); 1076dae52e14SToby Isaac PetscFunctionReturn(0); 1077dae52e14SToby Isaac } 1078dae52e14SToby Isaac 1079d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Gather(PetscPartitioner part, PetscViewer viewer) 1080dae52e14SToby Isaac { 1081dae52e14SToby Isaac PetscBool iascii; 1082dae52e14SToby Isaac PetscErrorCode ierr; 1083dae52e14SToby Isaac 1084dae52e14SToby Isaac PetscFunctionBegin; 1085dae52e14SToby Isaac PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 1086dae52e14SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1087dae52e14SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 1088dae52e14SToby Isaac if (iascii) {ierr = PetscPartitionerView_Gather_Ascii(part, viewer);CHKERRQ(ierr);} 1089dae52e14SToby Isaac PetscFunctionReturn(0); 1090dae52e14SToby Isaac } 1091dae52e14SToby Isaac 1092d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Gather(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 1093dae52e14SToby Isaac { 1094dae52e14SToby Isaac PetscInt np; 1095dae52e14SToby Isaac PetscErrorCode ierr; 1096dae52e14SToby Isaac 1097dae52e14SToby Isaac PetscFunctionBegin; 1098dae52e14SToby Isaac ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 1099dae52e14SToby Isaac ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr); 1100dae52e14SToby Isaac ierr = PetscSectionSetDof(partSection,0,numVertices);CHKERRQ(ierr); 1101dae52e14SToby Isaac for (np = 1; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, 0);CHKERRQ(ierr);} 1102dae52e14SToby Isaac ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 1103dae52e14SToby Isaac PetscFunctionReturn(0); 1104dae52e14SToby Isaac } 1105dae52e14SToby Isaac 1106d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Gather(PetscPartitioner part) 1107dae52e14SToby Isaac { 1108dae52e14SToby Isaac PetscFunctionBegin; 1109dae52e14SToby Isaac part->ops->view = PetscPartitionerView_Gather; 1110dae52e14SToby Isaac part->ops->destroy = PetscPartitionerDestroy_Gather; 1111dae52e14SToby Isaac part->ops->partition = PetscPartitionerPartition_Gather; 1112dae52e14SToby Isaac PetscFunctionReturn(0); 1113dae52e14SToby Isaac } 1114dae52e14SToby Isaac 1115dae52e14SToby Isaac /*MC 1116dae52e14SToby Isaac PETSCPARTITIONERGATHER = "gather" - A PetscPartitioner object 1117dae52e14SToby Isaac 1118dae52e14SToby Isaac Level: intermediate 1119dae52e14SToby Isaac 1120dae52e14SToby Isaac .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 1121dae52e14SToby Isaac M*/ 1122dae52e14SToby Isaac 1123dae52e14SToby Isaac PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Gather(PetscPartitioner part) 1124dae52e14SToby Isaac { 1125dae52e14SToby Isaac PetscPartitioner_Gather *p; 1126dae52e14SToby Isaac PetscErrorCode ierr; 1127dae52e14SToby Isaac 1128dae52e14SToby Isaac PetscFunctionBegin; 1129dae52e14SToby Isaac PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 1130dae52e14SToby Isaac ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 1131dae52e14SToby Isaac part->data = p; 1132dae52e14SToby Isaac 1133dae52e14SToby Isaac ierr = PetscPartitionerInitialize_Gather(part);CHKERRQ(ierr); 1134dae52e14SToby Isaac PetscFunctionReturn(0); 1135dae52e14SToby Isaac } 1136dae52e14SToby Isaac 1137dae52e14SToby Isaac 1138d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part) 113977623264SMatthew G. Knepley { 114077623264SMatthew G. Knepley PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data; 114177623264SMatthew G. Knepley PetscErrorCode ierr; 114277623264SMatthew G. Knepley 114377623264SMatthew G. Knepley PetscFunctionBegin; 114477623264SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 114577623264SMatthew G. Knepley PetscFunctionReturn(0); 114677623264SMatthew G. Knepley } 114777623264SMatthew G. Knepley 1148d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Chaco_Ascii(PetscPartitioner part, PetscViewer viewer) 114977623264SMatthew G. Knepley { 115077623264SMatthew G. Knepley PetscViewerFormat format; 115177623264SMatthew G. Knepley PetscErrorCode ierr; 115277623264SMatthew G. Knepley 115377623264SMatthew G. Knepley PetscFunctionBegin; 115477623264SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 115577623264SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "Chaco Graph Partitioner:\n");CHKERRQ(ierr); 115677623264SMatthew G. Knepley PetscFunctionReturn(0); 115777623264SMatthew G. Knepley } 115877623264SMatthew G. Knepley 1159d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer) 116077623264SMatthew G. Knepley { 116177623264SMatthew G. Knepley PetscBool iascii; 116277623264SMatthew G. Knepley PetscErrorCode ierr; 116377623264SMatthew G. Knepley 116477623264SMatthew G. Knepley PetscFunctionBegin; 116577623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 116677623264SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 116777623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 116877623264SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_Chaco_Ascii(part, viewer);CHKERRQ(ierr);} 116977623264SMatthew G. Knepley PetscFunctionReturn(0); 117077623264SMatthew G. Knepley } 117177623264SMatthew G. Knepley 117270034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO) 117370034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H) 117470034214SMatthew G. Knepley #include <unistd.h> 117570034214SMatthew G. Knepley #endif 117611d1e910SBarry Smith #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT) 117711d1e910SBarry Smith #include <chaco.h> 117811d1e910SBarry Smith #else 117911d1e910SBarry Smith /* Older versions of Chaco do not have an include file */ 118070034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts, 118170034214SMatthew G. Knepley float *ewgts, float *x, float *y, float *z, char *outassignname, 118270034214SMatthew G. Knepley char *outfilename, short *assignment, int architecture, int ndims_tot, 118370034214SMatthew G. Knepley int mesh_dims[3], double *goal, int global_method, int local_method, 118470034214SMatthew G. Knepley int rqi_flag, int vmax, int ndims, double eigtol, long seed); 118511d1e910SBarry Smith #endif 118670034214SMatthew G. Knepley extern int FREE_GRAPH; 118777623264SMatthew G. Knepley #endif 118870034214SMatthew G. Knepley 1189d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 119070034214SMatthew G. Knepley { 119177623264SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO) 119270034214SMatthew G. Knepley enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3}; 119370034214SMatthew G. Knepley MPI_Comm comm; 119470034214SMatthew G. Knepley int nvtxs = numVertices; /* number of vertices in full graph */ 119570034214SMatthew G. Knepley int *vwgts = NULL; /* weights for all vertices */ 119670034214SMatthew G. Knepley float *ewgts = NULL; /* weights for all edges */ 119770034214SMatthew G. Knepley float *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */ 119870034214SMatthew G. Knepley char *outassignname = NULL; /* name of assignment output file */ 119970034214SMatthew G. Knepley char *outfilename = NULL; /* output file name */ 120070034214SMatthew G. Knepley int architecture = 1; /* 0 => hypercube, d => d-dimensional mesh */ 120170034214SMatthew G. Knepley int ndims_tot = 0; /* total number of cube dimensions to divide */ 120270034214SMatthew G. Knepley int mesh_dims[3]; /* dimensions of mesh of processors */ 120370034214SMatthew G. Knepley double *goal = NULL; /* desired set sizes for each set */ 120470034214SMatthew G. Knepley int global_method = 1; /* global partitioning algorithm */ 120570034214SMatthew G. Knepley int local_method = 1; /* local partitioning algorithm */ 120670034214SMatthew G. Knepley int rqi_flag = 0; /* should I use RQI/Symmlq eigensolver? */ 120770034214SMatthew G. Knepley int vmax = 200; /* how many vertices to coarsen down to? */ 120870034214SMatthew G. Knepley int ndims = 1; /* number of eigenvectors (2^d sets) */ 120970034214SMatthew G. Knepley double eigtol = 0.001; /* tolerance on eigenvectors */ 121070034214SMatthew G. Knepley long seed = 123636512; /* for random graph mutations */ 121111d1e910SBarry Smith #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT) 121211d1e910SBarry Smith int *assignment; /* Output partition */ 121311d1e910SBarry Smith #else 121470034214SMatthew G. Knepley short int *assignment; /* Output partition */ 121511d1e910SBarry Smith #endif 121670034214SMatthew G. Knepley int fd_stdout, fd_pipe[2]; 121770034214SMatthew G. Knepley PetscInt *points; 121870034214SMatthew G. Knepley int i, v, p; 121970034214SMatthew G. Knepley PetscErrorCode ierr; 122070034214SMatthew G. Knepley 122170034214SMatthew G. Knepley PetscFunctionBegin; 122270034214SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 122370034214SMatthew G. Knepley if (!numVertices) { 122477623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 122577623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 122670034214SMatthew G. Knepley ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 122770034214SMatthew G. Knepley PetscFunctionReturn(0); 122870034214SMatthew G. Knepley } 122970034214SMatthew G. Knepley FREE_GRAPH = 0; /* Do not let Chaco free my memory */ 123070034214SMatthew G. Knepley for (i = 0; i < start[numVertices]; ++i) ++adjacency[i]; 123170034214SMatthew G. Knepley 123270034214SMatthew G. Knepley if (global_method == INERTIAL_METHOD) { 123370034214SMatthew G. Knepley /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */ 123470034214SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported"); 123570034214SMatthew G. Knepley } 123677623264SMatthew G. Knepley mesh_dims[0] = nparts; 123770034214SMatthew G. Knepley mesh_dims[1] = 1; 123870034214SMatthew G. Knepley mesh_dims[2] = 1; 123970034214SMatthew G. Knepley ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr); 124070034214SMatthew G. Knepley /* Chaco outputs to stdout. We redirect this to a buffer. */ 124170034214SMatthew G. Knepley /* TODO: check error codes for UNIX calls */ 124270034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H) 124370034214SMatthew G. Knepley { 124470034214SMatthew G. Knepley int piperet; 124570034214SMatthew G. Knepley piperet = pipe(fd_pipe); 124670034214SMatthew G. Knepley if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe"); 124770034214SMatthew G. Knepley fd_stdout = dup(1); 124870034214SMatthew G. Knepley close(1); 124970034214SMatthew G. Knepley dup2(fd_pipe[1], 1); 125070034214SMatthew G. Knepley } 125170034214SMatthew G. Knepley #endif 125270034214SMatthew G. Knepley ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename, 125370034214SMatthew G. Knepley assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag, 125470034214SMatthew G. Knepley vmax, ndims, eigtol, seed); 125570034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H) 125670034214SMatthew G. Knepley { 125770034214SMatthew G. Knepley char msgLog[10000]; 125870034214SMatthew G. Knepley int count; 125970034214SMatthew G. Knepley 126070034214SMatthew G. Knepley fflush(stdout); 126170034214SMatthew G. Knepley count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char)); 126270034214SMatthew G. Knepley if (count < 0) count = 0; 126370034214SMatthew G. Knepley msgLog[count] = 0; 126470034214SMatthew G. Knepley close(1); 126570034214SMatthew G. Knepley dup2(fd_stdout, 1); 126670034214SMatthew G. Knepley close(fd_stdout); 126770034214SMatthew G. Knepley close(fd_pipe[0]); 126870034214SMatthew G. Knepley close(fd_pipe[1]); 126970034214SMatthew G. Knepley if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog); 127070034214SMatthew G. Knepley } 127170034214SMatthew G. Knepley #endif 127270034214SMatthew G. Knepley /* Convert to PetscSection+IS */ 127377623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 127470034214SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) { 127577623264SMatthew G. Knepley ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr); 127670034214SMatthew G. Knepley } 127777623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 127870034214SMatthew G. Knepley ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr); 127977623264SMatthew G. Knepley for (p = 0, i = 0; p < nparts; ++p) { 128070034214SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) { 128170034214SMatthew G. Knepley if (assignment[v] == p) points[i++] = v; 128270034214SMatthew G. Knepley } 128370034214SMatthew G. Knepley } 128470034214SMatthew G. Knepley if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs); 128570034214SMatthew G. Knepley ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 128670034214SMatthew G. Knepley if (global_method == INERTIAL_METHOD) { 128770034214SMatthew G. Knepley /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */ 128870034214SMatthew G. Knepley } 128970034214SMatthew G. Knepley ierr = PetscFree(assignment);CHKERRQ(ierr); 129070034214SMatthew G. Knepley for (i = 0; i < start[numVertices]; ++i) --adjacency[i]; 129170034214SMatthew G. Knepley PetscFunctionReturn(0); 129277623264SMatthew G. Knepley #else 129377623264SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco."); 129470034214SMatthew G. Knepley #endif 129577623264SMatthew G. Knepley } 129677623264SMatthew G. Knepley 1297d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part) 129877623264SMatthew G. Knepley { 129977623264SMatthew G. Knepley PetscFunctionBegin; 130077623264SMatthew G. Knepley part->ops->view = PetscPartitionerView_Chaco; 130177623264SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_Chaco; 130277623264SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_Chaco; 130377623264SMatthew G. Knepley PetscFunctionReturn(0); 130477623264SMatthew G. Knepley } 130577623264SMatthew G. Knepley 130677623264SMatthew G. Knepley /*MC 130777623264SMatthew G. Knepley PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library 130877623264SMatthew G. Knepley 130977623264SMatthew G. Knepley Level: intermediate 131077623264SMatthew G. Knepley 131177623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 131277623264SMatthew G. Knepley M*/ 131377623264SMatthew G. Knepley 131477623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part) 131577623264SMatthew G. Knepley { 131677623264SMatthew G. Knepley PetscPartitioner_Chaco *p; 131777623264SMatthew G. Knepley PetscErrorCode ierr; 131877623264SMatthew G. Knepley 131977623264SMatthew G. Knepley PetscFunctionBegin; 132077623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 132177623264SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 132277623264SMatthew G. Knepley part->data = p; 132377623264SMatthew G. Knepley 132477623264SMatthew G. Knepley ierr = PetscPartitionerInitialize_Chaco(part);CHKERRQ(ierr); 132577623264SMatthew G. Knepley ierr = PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionercite);CHKERRQ(ierr); 132677623264SMatthew G. Knepley PetscFunctionReturn(0); 132777623264SMatthew G. Knepley } 132877623264SMatthew G. Knepley 1329d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part) 133077623264SMatthew G. Knepley { 133177623264SMatthew G. Knepley PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *) part->data; 133277623264SMatthew G. Knepley PetscErrorCode ierr; 133377623264SMatthew G. Knepley 133477623264SMatthew G. Knepley PetscFunctionBegin; 133577623264SMatthew G. Knepley ierr = PetscFree(p);CHKERRQ(ierr); 133677623264SMatthew G. Knepley PetscFunctionReturn(0); 133777623264SMatthew G. Knepley } 133877623264SMatthew G. Knepley 1339d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_ParMetis_Ascii(PetscPartitioner part, PetscViewer viewer) 134077623264SMatthew G. Knepley { 134177623264SMatthew G. Knepley PetscViewerFormat format; 134277623264SMatthew G. Knepley PetscErrorCode ierr; 134377623264SMatthew G. Knepley 134477623264SMatthew G. Knepley PetscFunctionBegin; 134577623264SMatthew G. Knepley ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 134677623264SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "ParMetis Graph Partitioner:\n");CHKERRQ(ierr); 134777623264SMatthew G. Knepley PetscFunctionReturn(0); 134877623264SMatthew G. Knepley } 134977623264SMatthew G. Knepley 1350d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer) 135177623264SMatthew G. Knepley { 135277623264SMatthew G. Knepley PetscBool iascii; 135377623264SMatthew G. Knepley PetscErrorCode ierr; 135477623264SMatthew G. Knepley 135577623264SMatthew G. Knepley PetscFunctionBegin; 135677623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 135777623264SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 135877623264SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 135977623264SMatthew G. Knepley if (iascii) {ierr = PetscPartitionerView_ParMetis_Ascii(part, viewer);CHKERRQ(ierr);} 136077623264SMatthew G. Knepley PetscFunctionReturn(0); 136177623264SMatthew G. Knepley } 136270034214SMatthew G. Knepley 136370034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS) 136470034214SMatthew G. Knepley #include <parmetis.h> 136577623264SMatthew G. Knepley #endif 136670034214SMatthew G. Knepley 1367d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition) 136870034214SMatthew G. Knepley { 136977623264SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS) 137070034214SMatthew G. Knepley MPI_Comm comm; 1371deea36a5SMatthew G. Knepley PetscSection section; 137270034214SMatthew G. Knepley PetscInt nvtxs = numVertices; /* The number of vertices in full graph */ 137370034214SMatthew G. Knepley PetscInt *vtxdist; /* Distribution of vertices across processes */ 137470034214SMatthew G. Knepley PetscInt *xadj = start; /* Start of edge list for each vertex */ 137570034214SMatthew G. Knepley PetscInt *adjncy = adjacency; /* Edge lists for all vertices */ 137670034214SMatthew G. Knepley PetscInt *vwgt = NULL; /* Vertex weights */ 137770034214SMatthew G. Knepley PetscInt *adjwgt = NULL; /* Edge weights */ 137870034214SMatthew G. Knepley PetscInt wgtflag = 0; /* Indicates which weights are present */ 137970034214SMatthew G. Knepley PetscInt numflag = 0; /* Indicates initial offset (0 or 1) */ 138070034214SMatthew G. Knepley PetscInt ncon = 1; /* The number of weights per vertex */ 138170034214SMatthew G. Knepley PetscReal *tpwgts; /* The fraction of vertex weights assigned to each partition */ 138270034214SMatthew G. Knepley PetscReal *ubvec; /* The balance intolerance for vertex weights */ 138370034214SMatthew G. Knepley PetscInt options[5]; /* Options */ 138470034214SMatthew G. Knepley /* Outputs */ 138570034214SMatthew G. Knepley PetscInt edgeCut; /* The number of edges cut by the partition */ 138670034214SMatthew G. Knepley PetscInt *assignment, *points; 138777623264SMatthew G. Knepley PetscMPIInt rank, p, v, i; 138870034214SMatthew G. Knepley PetscErrorCode ierr; 138970034214SMatthew G. Knepley 139070034214SMatthew G. Knepley PetscFunctionBegin; 139177623264SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr); 139270034214SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 139370034214SMatthew G. Knepley options[0] = 0; /* Use all defaults */ 139470034214SMatthew G. Knepley /* Calculate vertex distribution */ 1395cd0de0f2SShri ierr = PetscMalloc5(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment,nvtxs,&vwgt);CHKERRQ(ierr); 139670034214SMatthew G. Knepley vtxdist[0] = 0; 139770034214SMatthew G. Knepley ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr); 139870034214SMatthew G. Knepley for (p = 2; p <= nparts; ++p) { 139970034214SMatthew G. Knepley vtxdist[p] += vtxdist[p-1]; 140070034214SMatthew G. Knepley } 140170034214SMatthew G. Knepley /* Calculate weights */ 140270034214SMatthew G. Knepley for (p = 0; p < nparts; ++p) { 140370034214SMatthew G. Knepley tpwgts[p] = 1.0/nparts; 140470034214SMatthew G. Knepley } 140570034214SMatthew G. Knepley ubvec[0] = 1.05; 1406deea36a5SMatthew G. Knepley /* Weight cells by dofs on cell by default */ 1407deea36a5SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 1408deea36a5SMatthew G. Knepley if (section) { 1409deea36a5SMatthew G. Knepley PetscInt cStart, cEnd, dof; 141070034214SMatthew G. Knepley 1411deea36a5SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 1412deea36a5SMatthew G. Knepley for (v = cStart; v < cEnd; ++v) { 1413deea36a5SMatthew G. Knepley ierr = PetscSectionGetDof(section, v, &dof);CHKERRQ(ierr); 1414deea36a5SMatthew G. Knepley vwgt[v-cStart] = PetscMax(dof, 1); 1415deea36a5SMatthew G. Knepley } 1416deea36a5SMatthew G. Knepley } else { 1417deea36a5SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) vwgt[v] = 1; 1418cd0de0f2SShri } 1419cd0de0f2SShri 142070034214SMatthew G. Knepley if (nparts == 1) { 14219fc93327SToby Isaac ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));CHKERRQ(ierr); 142270034214SMatthew G. Knepley } else { 142370034214SMatthew G. Knepley if (vtxdist[1] == vtxdist[nparts]) { 142470034214SMatthew G. Knepley if (!rank) { 142570034214SMatthew G. Knepley PetscStackPush("METIS_PartGraphKway"); 142670034214SMatthew G. Knepley ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment); 142770034214SMatthew G. Knepley PetscStackPop; 142870034214SMatthew G. Knepley if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()"); 142970034214SMatthew G. Knepley } 143070034214SMatthew G. Knepley } else { 143170034214SMatthew G. Knepley PetscStackPush("ParMETIS_V3_PartKway"); 143270034214SMatthew G. Knepley ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm); 143370034214SMatthew G. Knepley PetscStackPop; 1434c717d290SMatthew G. Knepley if (ierr != METIS_OK) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error %d in ParMETIS_V3_PartKway()", ierr); 143570034214SMatthew G. Knepley } 143670034214SMatthew G. Knepley } 143770034214SMatthew G. Knepley /* Convert to PetscSection+IS */ 143877623264SMatthew G. Knepley ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr); 143977623264SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);} 144077623264SMatthew G. Knepley ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr); 144170034214SMatthew G. Knepley ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr); 144277623264SMatthew G. Knepley for (p = 0, i = 0; p < nparts; ++p) { 144370034214SMatthew G. Knepley for (v = 0; v < nvtxs; ++v) { 144470034214SMatthew G. Knepley if (assignment[v] == p) points[i++] = v; 144570034214SMatthew G. Knepley } 144670034214SMatthew G. Knepley } 144770034214SMatthew G. Knepley if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs); 144870034214SMatthew G. Knepley ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr); 1449cd0de0f2SShri ierr = PetscFree5(vtxdist,tpwgts,ubvec,assignment,vwgt);CHKERRQ(ierr); 14509b80ac48SMatthew G. Knepley PetscFunctionReturn(0); 145170034214SMatthew G. Knepley #else 145277623264SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis."); 145370034214SMatthew G. Knepley #endif 145470034214SMatthew G. Knepley } 145570034214SMatthew G. Knepley 1456d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part) 145777623264SMatthew G. Knepley { 145877623264SMatthew G. Knepley PetscFunctionBegin; 145977623264SMatthew G. Knepley part->ops->view = PetscPartitionerView_ParMetis; 146077623264SMatthew G. Knepley part->ops->destroy = PetscPartitionerDestroy_ParMetis; 146177623264SMatthew G. Knepley part->ops->partition = PetscPartitionerPartition_ParMetis; 146277623264SMatthew G. Knepley PetscFunctionReturn(0); 146377623264SMatthew G. Knepley } 146477623264SMatthew G. Knepley 146577623264SMatthew G. Knepley /*MC 146677623264SMatthew G. Knepley PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMetis library 146777623264SMatthew G. Knepley 146877623264SMatthew G. Knepley Level: intermediate 146977623264SMatthew G. Knepley 147077623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 147177623264SMatthew G. Knepley M*/ 147277623264SMatthew G. Knepley 147377623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part) 147477623264SMatthew G. Knepley { 147577623264SMatthew G. Knepley PetscPartitioner_ParMetis *p; 147677623264SMatthew G. Knepley PetscErrorCode ierr; 147777623264SMatthew G. Knepley 147877623264SMatthew G. Knepley PetscFunctionBegin; 147977623264SMatthew G. Knepley PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 148077623264SMatthew G. Knepley ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 148177623264SMatthew G. Knepley part->data = p; 148277623264SMatthew G. Knepley 148377623264SMatthew G. Knepley ierr = PetscPartitionerInitialize_ParMetis(part);CHKERRQ(ierr); 148477623264SMatthew G. Knepley ierr = PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionercite);CHKERRQ(ierr); 148570034214SMatthew G. Knepley PetscFunctionReturn(0); 148670034214SMatthew G. Knepley } 148770034214SMatthew G. Knepley 14885680f57bSMatthew G. Knepley /*@ 14895680f57bSMatthew G. Knepley DMPlexGetPartitioner - Get the mesh partitioner 14905680f57bSMatthew G. Knepley 14915680f57bSMatthew G. Knepley Not collective 14925680f57bSMatthew G. Knepley 14935680f57bSMatthew G. Knepley Input Parameter: 14945680f57bSMatthew G. Knepley . dm - The DM 14955680f57bSMatthew G. Knepley 14965680f57bSMatthew G. Knepley Output Parameter: 14975680f57bSMatthew G. Knepley . part - The PetscPartitioner 14985680f57bSMatthew G. Knepley 14995680f57bSMatthew G. Knepley Level: developer 15005680f57bSMatthew G. Knepley 150198599a47SLawrence Mitchell Note: This gets a borrowed reference, so the user should not destroy this PetscPartitioner. 150298599a47SLawrence Mitchell 150398599a47SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexSetPartitioner(), PetscPartitionerCreate() 15045680f57bSMatthew G. Knepley @*/ 15055680f57bSMatthew G. Knepley PetscErrorCode DMPlexGetPartitioner(DM dm, PetscPartitioner *part) 15065680f57bSMatthew G. Knepley { 15075680f57bSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 15084f3833eaSMatthew G. Knepley PetscErrorCode ierr; 15095680f57bSMatthew G. Knepley 15105680f57bSMatthew G. Knepley PetscFunctionBegin; 15115680f57bSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15125680f57bSMatthew G. Knepley PetscValidPointer(part, 2); 15134f3833eaSMatthew G. Knepley ierr = PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);CHKERRQ(ierr); 15145680f57bSMatthew G. Knepley *part = mesh->partitioner; 15155680f57bSMatthew G. Knepley PetscFunctionReturn(0); 15165680f57bSMatthew G. Knepley } 15175680f57bSMatthew G. Knepley 151871bb2955SLawrence Mitchell /*@ 151971bb2955SLawrence Mitchell DMPlexSetPartitioner - Set the mesh partitioner 152071bb2955SLawrence Mitchell 152171bb2955SLawrence Mitchell logically collective on dm and part 152271bb2955SLawrence Mitchell 152371bb2955SLawrence Mitchell Input Parameters: 152471bb2955SLawrence Mitchell + dm - The DM 152571bb2955SLawrence Mitchell - part - The partitioner 152671bb2955SLawrence Mitchell 152771bb2955SLawrence Mitchell Level: developer 152871bb2955SLawrence Mitchell 152971bb2955SLawrence Mitchell Note: Any existing PetscPartitioner will be destroyed. 153071bb2955SLawrence Mitchell 153171bb2955SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexGetPartitioner(), PetscPartitionerCreate() 153271bb2955SLawrence Mitchell @*/ 153371bb2955SLawrence Mitchell PetscErrorCode DMPlexSetPartitioner(DM dm, PetscPartitioner part) 153471bb2955SLawrence Mitchell { 153571bb2955SLawrence Mitchell DM_Plex *mesh = (DM_Plex *) dm->data; 153671bb2955SLawrence Mitchell PetscErrorCode ierr; 153771bb2955SLawrence Mitchell 153871bb2955SLawrence Mitchell PetscFunctionBegin; 153971bb2955SLawrence Mitchell PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 154071bb2955SLawrence Mitchell PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 2); 154171bb2955SLawrence Mitchell ierr = PetscObjectReference((PetscObject)part);CHKERRQ(ierr); 154271bb2955SLawrence Mitchell ierr = PetscPartitionerDestroy(&mesh->partitioner);CHKERRQ(ierr); 154371bb2955SLawrence Mitchell mesh->partitioner = part; 154471bb2955SLawrence Mitchell PetscFunctionReturn(0); 154571bb2955SLawrence Mitchell } 154671bb2955SLawrence Mitchell 15476a5a2ffdSToby Isaac static PetscErrorCode DMPlexPartitionLabelClosure_Tree(DM dm, DMLabel label, PetscInt rank, PetscInt point, PetscBool up, PetscBool down) 1548270bba0cSToby Isaac { 1549270bba0cSToby Isaac PetscErrorCode ierr; 1550270bba0cSToby Isaac 1551270bba0cSToby Isaac PetscFunctionBegin; 15526a5a2ffdSToby Isaac if (up) { 15536a5a2ffdSToby Isaac PetscInt parent; 15546a5a2ffdSToby Isaac 1555270bba0cSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr); 15566a5a2ffdSToby Isaac if (parent != point) { 15576a5a2ffdSToby Isaac PetscInt closureSize, *closure = NULL, i; 15586a5a2ffdSToby Isaac 1559270bba0cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 1560270bba0cSToby Isaac for (i = 0; i < closureSize; i++) { 1561270bba0cSToby Isaac PetscInt cpoint = closure[2*i]; 1562270bba0cSToby Isaac 1563270bba0cSToby Isaac ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr); 15646a5a2ffdSToby Isaac ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 1565270bba0cSToby Isaac } 1566270bba0cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 15676a5a2ffdSToby Isaac } 15686a5a2ffdSToby Isaac } 15696a5a2ffdSToby Isaac if (down) { 15706a5a2ffdSToby Isaac PetscInt numChildren; 15716a5a2ffdSToby Isaac const PetscInt *children; 15726a5a2ffdSToby Isaac 15736a5a2ffdSToby Isaac ierr = DMPlexGetTreeChildren(dm,point,&numChildren,&children);CHKERRQ(ierr); 15746a5a2ffdSToby Isaac if (numChildren) { 15756a5a2ffdSToby Isaac PetscInt i; 15766a5a2ffdSToby Isaac 15776a5a2ffdSToby Isaac for (i = 0; i < numChildren; i++) { 15786a5a2ffdSToby Isaac PetscInt cpoint = children[i]; 15796a5a2ffdSToby Isaac 15806a5a2ffdSToby Isaac ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr); 15816a5a2ffdSToby Isaac ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr); 15826a5a2ffdSToby Isaac } 15836a5a2ffdSToby Isaac } 15846a5a2ffdSToby Isaac } 1585270bba0cSToby Isaac PetscFunctionReturn(0); 1586270bba0cSToby Isaac } 1587270bba0cSToby Isaac 15885abbe4feSMichael Lange /*@ 15895abbe4feSMichael Lange DMPlexPartitionLabelClosure - Add the closure of all points to the partition label 15905abbe4feSMichael Lange 15915abbe4feSMichael Lange Input Parameters: 15925abbe4feSMichael Lange + dm - The DM 15935abbe4feSMichael Lange - label - DMLabel assinging ranks to remote roots 15945abbe4feSMichael Lange 15955abbe4feSMichael Lange Level: developer 15965abbe4feSMichael Lange 15975abbe4feSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 15985abbe4feSMichael Lange @*/ 15995abbe4feSMichael Lange PetscErrorCode DMPlexPartitionLabelClosure(DM dm, DMLabel label) 16009ffc88e4SToby Isaac { 16015abbe4feSMichael Lange IS rankIS, pointIS; 16025abbe4feSMichael Lange const PetscInt *ranks, *points; 16035abbe4feSMichael Lange PetscInt numRanks, numPoints, r, p, c, closureSize; 16045abbe4feSMichael Lange PetscInt *closure = NULL; 16059ffc88e4SToby Isaac PetscErrorCode ierr; 16069ffc88e4SToby Isaac 16079ffc88e4SToby Isaac PetscFunctionBegin; 16085abbe4feSMichael Lange ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr); 16095abbe4feSMichael Lange ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 16105abbe4feSMichael Lange ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 16115abbe4feSMichael Lange for (r = 0; r < numRanks; ++r) { 16125abbe4feSMichael Lange const PetscInt rank = ranks[r]; 16139ffc88e4SToby Isaac 16145abbe4feSMichael Lange ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr); 16155abbe4feSMichael Lange ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 16165abbe4feSMichael Lange ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 16175abbe4feSMichael Lange for (p = 0; p < numPoints; ++p) { 16185abbe4feSMichael Lange ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 1619270bba0cSToby Isaac for (c = 0; c < closureSize*2; c += 2) { 1620270bba0cSToby Isaac ierr = DMLabelSetValue(label, closure[c], rank);CHKERRQ(ierr); 16216a5a2ffdSToby Isaac ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,closure[c],PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr); 1622270bba0cSToby Isaac } 16239ffc88e4SToby Isaac } 16245abbe4feSMichael Lange ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 16255abbe4feSMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 16269ffc88e4SToby Isaac } 16277de78196SMichael Lange if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);} 16285abbe4feSMichael Lange ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 16295abbe4feSMichael Lange ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 16309ffc88e4SToby Isaac PetscFunctionReturn(0); 16319ffc88e4SToby Isaac } 16329ffc88e4SToby Isaac 163324d039d7SMichael Lange /*@ 163424d039d7SMichael Lange DMPlexPartitionLabelAdjacency - Add one level of adjacent points to the partition label 163524d039d7SMichael Lange 163624d039d7SMichael Lange Input Parameters: 163724d039d7SMichael Lange + dm - The DM 163824d039d7SMichael Lange - label - DMLabel assinging ranks to remote roots 163924d039d7SMichael Lange 164024d039d7SMichael Lange Level: developer 164124d039d7SMichael Lange 164224d039d7SMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 164324d039d7SMichael Lange @*/ 164424d039d7SMichael Lange PetscErrorCode DMPlexPartitionLabelAdjacency(DM dm, DMLabel label) 164570034214SMatthew G. Knepley { 164624d039d7SMichael Lange IS rankIS, pointIS; 164724d039d7SMichael Lange const PetscInt *ranks, *points; 164824d039d7SMichael Lange PetscInt numRanks, numPoints, r, p, a, adjSize; 164924d039d7SMichael Lange PetscInt *adj = NULL; 165070034214SMatthew G. Knepley PetscErrorCode ierr; 165170034214SMatthew G. Knepley 165270034214SMatthew G. Knepley PetscFunctionBegin; 165324d039d7SMichael Lange ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr); 165424d039d7SMichael Lange ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 165524d039d7SMichael Lange ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 165624d039d7SMichael Lange for (r = 0; r < numRanks; ++r) { 165724d039d7SMichael Lange const PetscInt rank = ranks[r]; 165870034214SMatthew G. Knepley 165924d039d7SMichael Lange ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr); 166024d039d7SMichael Lange ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 166124d039d7SMichael Lange ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 166270034214SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 166324d039d7SMichael Lange adjSize = PETSC_DETERMINE; 166424d039d7SMichael Lange ierr = DMPlexGetAdjacency(dm, points[p], &adjSize, &adj);CHKERRQ(ierr); 166524d039d7SMichael Lange for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(label, adj[a], rank);CHKERRQ(ierr);} 166670034214SMatthew G. Knepley } 166724d039d7SMichael Lange ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 166824d039d7SMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 166970034214SMatthew G. Knepley } 167024d039d7SMichael Lange ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 167124d039d7SMichael Lange ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 167224d039d7SMichael Lange ierr = PetscFree(adj);CHKERRQ(ierr); 167324d039d7SMichael Lange PetscFunctionReturn(0); 167470034214SMatthew G. Knepley } 167570034214SMatthew G. Knepley 1676be200f8dSMichael Lange /*@ 1677be200f8dSMichael Lange DMPlexPartitionLabelPropagate - Propagate points in a partition label over the point SF 1678be200f8dSMichael Lange 1679be200f8dSMichael Lange Input Parameters: 1680be200f8dSMichael Lange + dm - The DM 1681be200f8dSMichael Lange - label - DMLabel assinging ranks to remote roots 1682be200f8dSMichael Lange 1683be200f8dSMichael Lange Level: developer 1684be200f8dSMichael Lange 1685be200f8dSMichael Lange Note: This is required when generating multi-level overlaps to capture 1686be200f8dSMichael Lange overlap points from non-neighbouring partitions. 1687be200f8dSMichael Lange 1688be200f8dSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 1689be200f8dSMichael Lange @*/ 1690be200f8dSMichael Lange PetscErrorCode DMPlexPartitionLabelPropagate(DM dm, DMLabel label) 1691be200f8dSMichael Lange { 1692be200f8dSMichael Lange MPI_Comm comm; 1693be200f8dSMichael Lange PetscMPIInt rank; 1694be200f8dSMichael Lange PetscSF sfPoint; 16955d04f6ebSMichael Lange DMLabel lblRoots, lblLeaves; 1696be200f8dSMichael Lange IS rankIS, pointIS; 1697be200f8dSMichael Lange const PetscInt *ranks; 1698be200f8dSMichael Lange PetscInt numRanks, r; 1699be200f8dSMichael Lange PetscErrorCode ierr; 1700be200f8dSMichael Lange 1701be200f8dSMichael Lange PetscFunctionBegin; 1702be200f8dSMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 1703be200f8dSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1704be200f8dSMichael Lange ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 17055d04f6ebSMichael Lange /* Pull point contributions from remote leaves into local roots */ 17065d04f6ebSMichael Lange ierr = DMLabelGather(label, sfPoint, &lblLeaves);CHKERRQ(ierr); 17075d04f6ebSMichael Lange ierr = DMLabelGetValueIS(lblLeaves, &rankIS);CHKERRQ(ierr); 17085d04f6ebSMichael Lange ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 17095d04f6ebSMichael Lange ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 17105d04f6ebSMichael Lange for (r = 0; r < numRanks; ++r) { 17115d04f6ebSMichael Lange const PetscInt remoteRank = ranks[r]; 17125d04f6ebSMichael Lange if (remoteRank == rank) continue; 17135d04f6ebSMichael Lange ierr = DMLabelGetStratumIS(lblLeaves, remoteRank, &pointIS);CHKERRQ(ierr); 17145d04f6ebSMichael Lange ierr = DMLabelInsertIS(label, pointIS, remoteRank);CHKERRQ(ierr); 17155d04f6ebSMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 17165d04f6ebSMichael Lange } 17175d04f6ebSMichael Lange ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 17185d04f6ebSMichael Lange ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 17195d04f6ebSMichael Lange ierr = DMLabelDestroy(&lblLeaves);CHKERRQ(ierr); 1720be200f8dSMichael Lange /* Push point contributions from roots into remote leaves */ 1721be200f8dSMichael Lange ierr = DMLabelDistribute(label, sfPoint, &lblRoots);CHKERRQ(ierr); 1722be200f8dSMichael Lange ierr = DMLabelGetValueIS(lblRoots, &rankIS);CHKERRQ(ierr); 1723be200f8dSMichael Lange ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr); 1724be200f8dSMichael Lange ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr); 1725be200f8dSMichael Lange for (r = 0; r < numRanks; ++r) { 1726be200f8dSMichael Lange const PetscInt remoteRank = ranks[r]; 1727be200f8dSMichael Lange if (remoteRank == rank) continue; 1728be200f8dSMichael Lange ierr = DMLabelGetStratumIS(lblRoots, remoteRank, &pointIS);CHKERRQ(ierr); 1729be200f8dSMichael Lange ierr = DMLabelInsertIS(label, pointIS, remoteRank);CHKERRQ(ierr); 1730be200f8dSMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 1731be200f8dSMichael Lange } 1732be200f8dSMichael Lange ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr); 1733be200f8dSMichael Lange ierr = ISDestroy(&rankIS);CHKERRQ(ierr); 1734be200f8dSMichael Lange ierr = DMLabelDestroy(&lblRoots);CHKERRQ(ierr); 1735be200f8dSMichael Lange PetscFunctionReturn(0); 1736be200f8dSMichael Lange } 1737be200f8dSMichael Lange 17381fd9873aSMichael Lange /*@ 17391fd9873aSMichael Lange DMPlexPartitionLabelInvert - Create a partition label of remote roots from a local root label 17401fd9873aSMichael Lange 17411fd9873aSMichael Lange Input Parameters: 17421fd9873aSMichael Lange + dm - The DM 17431fd9873aSMichael Lange . rootLabel - DMLabel assinging ranks to local roots 17441fd9873aSMichael Lange . processSF - A star forest mapping into the local index on each remote rank 17451fd9873aSMichael Lange 17461fd9873aSMichael Lange Output Parameter: 17471fd9873aSMichael Lange - leafLabel - DMLabel assinging ranks to remote roots 17481fd9873aSMichael Lange 17491fd9873aSMichael Lange Note: The rootLabel defines a send pattern by mapping local points to remote target ranks. The 17501fd9873aSMichael Lange resulting leafLabel is a receiver mapping of remote roots to their parent rank. 17511fd9873aSMichael Lange 17521fd9873aSMichael Lange Level: developer 17531fd9873aSMichael Lange 17541fd9873aSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap 17551fd9873aSMichael Lange @*/ 17561fd9873aSMichael Lange PetscErrorCode DMPlexPartitionLabelInvert(DM dm, DMLabel rootLabel, PetscSF processSF, DMLabel leafLabel) 17571fd9873aSMichael Lange { 17581fd9873aSMichael Lange MPI_Comm comm; 17599852e123SBarry Smith PetscMPIInt rank, size; 17609852e123SBarry Smith PetscInt p, n, numNeighbors, ssize, l, nleaves; 17611fd9873aSMichael Lange PetscSF sfPoint; 17621fd9873aSMichael Lange PetscSFNode *rootPoints, *leafPoints; 17631fd9873aSMichael Lange PetscSection rootSection, leafSection; 17641fd9873aSMichael Lange const PetscSFNode *remote; 17651fd9873aSMichael Lange const PetscInt *local, *neighbors; 17661fd9873aSMichael Lange IS valueIS; 17671fd9873aSMichael Lange PetscErrorCode ierr; 17681fd9873aSMichael Lange 17691fd9873aSMichael Lange PetscFunctionBegin; 17701fd9873aSMichael Lange ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 17711fd9873aSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 17729852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 17731fd9873aSMichael Lange ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 17741fd9873aSMichael Lange 17751fd9873aSMichael Lange /* Convert to (point, rank) and use actual owners */ 17761fd9873aSMichael Lange ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 17779852e123SBarry Smith ierr = PetscSectionSetChart(rootSection, 0, size);CHKERRQ(ierr); 17781fd9873aSMichael Lange ierr = DMLabelGetValueIS(rootLabel, &valueIS);CHKERRQ(ierr); 17791fd9873aSMichael Lange ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr); 17801fd9873aSMichael Lange ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr); 17811fd9873aSMichael Lange for (n = 0; n < numNeighbors; ++n) { 17821fd9873aSMichael Lange PetscInt numPoints; 17831fd9873aSMichael Lange 17841fd9873aSMichael Lange ierr = DMLabelGetStratumSize(rootLabel, neighbors[n], &numPoints);CHKERRQ(ierr); 17851fd9873aSMichael Lange ierr = PetscSectionAddDof(rootSection, neighbors[n], numPoints);CHKERRQ(ierr); 17861fd9873aSMichael Lange } 17871fd9873aSMichael Lange ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 17889852e123SBarry Smith ierr = PetscSectionGetStorageSize(rootSection, &ssize);CHKERRQ(ierr); 17899852e123SBarry Smith ierr = PetscMalloc1(ssize, &rootPoints);CHKERRQ(ierr); 17901fd9873aSMichael Lange ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 17911fd9873aSMichael Lange for (n = 0; n < numNeighbors; ++n) { 17921fd9873aSMichael Lange IS pointIS; 17931fd9873aSMichael Lange const PetscInt *points; 17941fd9873aSMichael Lange PetscInt off, numPoints, p; 17951fd9873aSMichael Lange 17961fd9873aSMichael Lange ierr = PetscSectionGetOffset(rootSection, neighbors[n], &off);CHKERRQ(ierr); 17971fd9873aSMichael Lange ierr = DMLabelGetStratumIS(rootLabel, neighbors[n], &pointIS);CHKERRQ(ierr); 17981fd9873aSMichael Lange ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 17991fd9873aSMichael Lange ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 18001fd9873aSMichael Lange for (p = 0; p < numPoints; ++p) { 1801f8987ae8SMichael Lange if (local) {ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);} 1802f8987ae8SMichael Lange else {l = -1;} 18031fd9873aSMichael Lange if (l >= 0) {rootPoints[off+p] = remote[l];} 18041fd9873aSMichael Lange else {rootPoints[off+p].index = points[p]; rootPoints[off+p].rank = rank;} 18051fd9873aSMichael Lange } 18061fd9873aSMichael Lange ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 18071fd9873aSMichael Lange ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 18081fd9873aSMichael Lange } 18091fd9873aSMichael Lange ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr); 18101fd9873aSMichael Lange ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 18111fd9873aSMichael Lange /* Communicate overlap */ 18121fd9873aSMichael Lange ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); 18131fd9873aSMichael Lange ierr = DMPlexDistributeData(dm, processSF, rootSection, MPIU_2INT, rootPoints, leafSection, (void**) &leafPoints);CHKERRQ(ierr); 18141fd9873aSMichael Lange /* Filter remote contributions (ovLeafPoints) into the overlapSF */ 18159852e123SBarry Smith ierr = PetscSectionGetStorageSize(leafSection, &ssize);CHKERRQ(ierr); 18169852e123SBarry Smith for (p = 0; p < ssize; p++) { 18171fd9873aSMichael Lange ierr = DMLabelSetValue(leafLabel, leafPoints[p].index, leafPoints[p].rank);CHKERRQ(ierr); 18181fd9873aSMichael Lange } 18191fd9873aSMichael Lange ierr = PetscFree(rootPoints);CHKERRQ(ierr); 18201fd9873aSMichael Lange ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 18211fd9873aSMichael Lange ierr = PetscFree(leafPoints);CHKERRQ(ierr); 18221fd9873aSMichael Lange ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 18231fd9873aSMichael Lange PetscFunctionReturn(0); 18241fd9873aSMichael Lange } 18251fd9873aSMichael Lange 1826aa3148a8SMichael Lange /*@ 1827aa3148a8SMichael Lange DMPlexPartitionLabelCreateSF - Create a star forest from a label that assigns ranks to points 1828aa3148a8SMichael Lange 1829aa3148a8SMichael Lange Input Parameters: 1830aa3148a8SMichael Lange + dm - The DM 1831aa3148a8SMichael Lange . label - DMLabel assinging ranks to remote roots 1832aa3148a8SMichael Lange 1833aa3148a8SMichael Lange Output Parameter: 1834aa3148a8SMichael Lange - sf - The star forest communication context encapsulating the defined mapping 1835aa3148a8SMichael Lange 1836aa3148a8SMichael Lange Note: The incoming label is a receiver mapping of remote points to their parent rank. 1837aa3148a8SMichael Lange 1838aa3148a8SMichael Lange Level: developer 1839aa3148a8SMichael Lange 1840aa3148a8SMichael Lange .seealso: DMPlexDistribute(), DMPlexCreateOverlap 1841aa3148a8SMichael Lange @*/ 1842aa3148a8SMichael Lange PetscErrorCode DMPlexPartitionLabelCreateSF(DM dm, DMLabel label, PetscSF *sf) 1843aa3148a8SMichael Lange { 18449852e123SBarry Smith PetscMPIInt rank, size; 184543f7d02bSMichael Lange PetscInt n, numRemote, p, numPoints, pStart, pEnd, idx = 0; 1846aa3148a8SMichael Lange PetscSFNode *remotePoints; 184743f7d02bSMichael Lange IS remoteRootIS; 184843f7d02bSMichael Lange const PetscInt *remoteRoots; 1849aa3148a8SMichael Lange PetscErrorCode ierr; 1850aa3148a8SMichael Lange 1851aa3148a8SMichael Lange PetscFunctionBegin; 185243f7d02bSMichael Lange ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 18539852e123SBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr); 1854aa3148a8SMichael Lange 18559852e123SBarry Smith for (numRemote = 0, n = 0; n < size; ++n) { 1856aa3148a8SMichael Lange ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr); 1857aa3148a8SMichael Lange numRemote += numPoints; 1858aa3148a8SMichael Lange } 1859aa3148a8SMichael Lange ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr); 186043f7d02bSMichael Lange /* Put owned points first */ 186143f7d02bSMichael Lange ierr = DMLabelGetStratumSize(label, rank, &numPoints);CHKERRQ(ierr); 186243f7d02bSMichael Lange if (numPoints > 0) { 186343f7d02bSMichael Lange ierr = DMLabelGetStratumIS(label, rank, &remoteRootIS);CHKERRQ(ierr); 186443f7d02bSMichael Lange ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 186543f7d02bSMichael Lange for (p = 0; p < numPoints; p++) { 186643f7d02bSMichael Lange remotePoints[idx].index = remoteRoots[p]; 186743f7d02bSMichael Lange remotePoints[idx].rank = rank; 186843f7d02bSMichael Lange idx++; 186943f7d02bSMichael Lange } 187043f7d02bSMichael Lange ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 187143f7d02bSMichael Lange ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr); 187243f7d02bSMichael Lange } 187343f7d02bSMichael Lange /* Now add remote points */ 18749852e123SBarry Smith for (n = 0; n < size; ++n) { 1875aa3148a8SMichael Lange ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr); 187643f7d02bSMichael Lange if (numPoints <= 0 || n == rank) continue; 1877aa3148a8SMichael Lange ierr = DMLabelGetStratumIS(label, n, &remoteRootIS);CHKERRQ(ierr); 1878aa3148a8SMichael Lange ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 1879aa3148a8SMichael Lange for (p = 0; p < numPoints; p++) { 1880aa3148a8SMichael Lange remotePoints[idx].index = remoteRoots[p]; 1881aa3148a8SMichael Lange remotePoints[idx].rank = n; 1882aa3148a8SMichael Lange idx++; 1883aa3148a8SMichael Lange } 1884aa3148a8SMichael Lange ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr); 1885aa3148a8SMichael Lange ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr); 1886aa3148a8SMichael Lange } 1887aa3148a8SMichael Lange ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), sf);CHKERRQ(ierr); 1888aa3148a8SMichael Lange ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1889aa3148a8SMichael Lange ierr = PetscSFSetGraph(*sf, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 189070034214SMatthew G. Knepley PetscFunctionReturn(0); 189170034214SMatthew G. Knepley } 1892