1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2412e9a14SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For PetscFEInterpolate_Static() */ 3*012bc364SMatthew G. Knepley 4*012bc364SMatthew G. Knepley #include <petscdmplextransform.h> 575d3a19aSMatthew G. Knepley #include <petscsf.h> 675d3a19aSMatthew G. Knepley 7963fc26aSMatthew G. Knepley /*@ 8963fc26aSMatthew G. Knepley DMPlexCreateProcessSF - Create an SF which just has process connectivity 9963fc26aSMatthew G. Knepley 10d083f849SBarry Smith Collective on dm 11963fc26aSMatthew G. Knepley 12963fc26aSMatthew G. Knepley Input Parameters: 13963fc26aSMatthew G. Knepley + dm - The DM 14963fc26aSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity 15963fc26aSMatthew G. Knepley 16963fc26aSMatthew G. Knepley Output Parameters: 17963fc26aSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL 18963fc26aSMatthew G. Knepley - sfProcess - An SF encoding the process connectivity, or NULL 19963fc26aSMatthew G. Knepley 20963fc26aSMatthew G. Knepley Level: developer 21963fc26aSMatthew G. Knepley 22963fc26aSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateTwoSidedProcessSF() 23963fc26aSMatthew G. Knepley @*/ 24963fc26aSMatthew G. Knepley PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess) 2575d3a19aSMatthew G. Knepley { 2675d3a19aSMatthew G. Knepley PetscInt numRoots, numLeaves, l; 2775d3a19aSMatthew G. Knepley const PetscInt *localPoints; 2875d3a19aSMatthew G. Knepley const PetscSFNode *remotePoints; 2975d3a19aSMatthew G. Knepley PetscInt *localPointsNew; 3075d3a19aSMatthew G. Knepley PetscSFNode *remotePointsNew; 3175d3a19aSMatthew G. Knepley PetscInt *ranks, *ranksNew; 329852e123SBarry Smith PetscMPIInt size; 3375d3a19aSMatthew G. Knepley PetscErrorCode ierr; 3475d3a19aSMatthew G. Knepley 3575d3a19aSMatthew G. Knepley PetscFunctionBegin; 36963fc26aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 37963fc26aSMatthew G. Knepley PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 38963fc26aSMatthew G. Knepley if (processRanks) {PetscValidPointer(processRanks, 3);} 39963fc26aSMatthew G. Knepley if (sfProcess) {PetscValidPointer(sfProcess, 4);} 4055b25c41SPierre Jolivet ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRMPI(ierr); 4175d3a19aSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);CHKERRQ(ierr); 42785e854fSJed Brown ierr = PetscMalloc1(numLeaves, &ranks);CHKERRQ(ierr); 4375d3a19aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 4475d3a19aSMatthew G. Knepley ranks[l] = remotePoints[l].rank; 4575d3a19aSMatthew G. Knepley } 4675d3a19aSMatthew G. Knepley ierr = PetscSortRemoveDupsInt(&numLeaves, ranks);CHKERRQ(ierr); 47785e854fSJed Brown ierr = PetscMalloc1(numLeaves, &ranksNew);CHKERRQ(ierr); 48785e854fSJed Brown ierr = PetscMalloc1(numLeaves, &localPointsNew);CHKERRQ(ierr); 49785e854fSJed Brown ierr = PetscMalloc1(numLeaves, &remotePointsNew);CHKERRQ(ierr); 5075d3a19aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 5175d3a19aSMatthew G. Knepley ranksNew[l] = ranks[l]; 5275d3a19aSMatthew G. Knepley localPointsNew[l] = l; 5375d3a19aSMatthew G. Knepley remotePointsNew[l].index = 0; 5475d3a19aSMatthew G. Knepley remotePointsNew[l].rank = ranksNew[l]; 5575d3a19aSMatthew G. Knepley } 5675d3a19aSMatthew G. Knepley ierr = PetscFree(ranks);CHKERRQ(ierr); 57963fc26aSMatthew G. Knepley if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);} 58963fc26aSMatthew G. Knepley else {ierr = PetscFree(ranksNew);CHKERRQ(ierr);} 59963fc26aSMatthew G. Knepley if (sfProcess) { 6075d3a19aSMatthew G. Knepley ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr); 61963fc26aSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *sfProcess, "Process SF");CHKERRQ(ierr); 6275d3a19aSMatthew G. Knepley ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr); 639852e123SBarry Smith ierr = PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 64963fc26aSMatthew G. Knepley } 6575d3a19aSMatthew G. Knepley PetscFunctionReturn(0); 6675d3a19aSMatthew G. Knepley } 6775d3a19aSMatthew G. Knepley 682389894bSMatthew G. Knepley /*@ 692389894bSMatthew G. Knepley DMPlexCreateCoarsePointIS - Creates an IS covering the coarse DM chart with the fine points as data 702389894bSMatthew G. Knepley 71*012bc364SMatthew G. Knepley Collective on dm 72*012bc364SMatthew G. Knepley 732389894bSMatthew G. Knepley Input Parameter: 742389894bSMatthew G. Knepley . dm - The coarse DM 752389894bSMatthew G. Knepley 762389894bSMatthew G. Knepley Output Parameter: 772389894bSMatthew G. Knepley . fpointIS - The IS of all the fine points which exist in the original coarse mesh 782389894bSMatthew G. Knepley 792389894bSMatthew G. Knepley Level: developer 802389894bSMatthew G. Knepley 8197d8846cSMatthew Knepley .seealso: DMRefine(), DMPlexSetRefinementUniform(), DMPlexGetSubpointIS() 822389894bSMatthew G. Knepley @*/ 832389894bSMatthew G. Knepley PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS) 842389894bSMatthew G. Knepley { 85*012bc364SMatthew G. Knepley DMPlexTransform tr; 86412e9a14SMatthew G. Knepley PetscInt *fpoints; 87327c2912SStefano Zampini PetscInt pStart, pEnd, p, vStart, vEnd, v; 882389894bSMatthew G. Knepley PetscErrorCode ierr; 892389894bSMatthew G. Knepley 902389894bSMatthew G. Knepley PetscFunctionBegin; 912389894bSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 922389894bSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 93*012bc364SMatthew G. Knepley ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr);CHKERRQ(ierr); 94*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetUp(tr);CHKERRQ(ierr); 952389894bSMatthew G. Knepley ierr = PetscMalloc1(pEnd-pStart, &fpoints);CHKERRQ(ierr); 962389894bSMatthew G. Knepley for (p = 0; p < pEnd-pStart; ++p) fpoints[p] = -1; 97412e9a14SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 98*012bc364SMatthew G. Knepley PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */ 99327c2912SStefano Zampini 100*012bc364SMatthew G. Knepley ierr = DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew);CHKERRQ(ierr); 101412e9a14SMatthew G. Knepley fpoints[v-pStart] = vNew; 1022389894bSMatthew G. Knepley } 103*012bc364SMatthew G. Knepley ierr = DMPlexTransformDestroy(&tr);CHKERRQ(ierr); 1042389894bSMatthew G. Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd-pStart, fpoints, PETSC_OWN_POINTER, fpointIS);CHKERRQ(ierr); 1052389894bSMatthew G. Knepley PetscFunctionReturn(0); 1062389894bSMatthew G. Knepley } 1072389894bSMatthew G. Knepley 108*012bc364SMatthew G. Knepley /*@C 109*012bc364SMatthew G. Knepley DMPlexSetTransformType - Set the transform type for uniform refinement 110*012bc364SMatthew G. Knepley 111*012bc364SMatthew G. Knepley Input Parameters: 112*012bc364SMatthew G. Knepley + dm - The DM 113*012bc364SMatthew G. Knepley - type - The transform type for uniform refinement 114*012bc364SMatthew G. Knepley 115*012bc364SMatthew G. Knepley Level: developer 116*012bc364SMatthew G. Knepley 117*012bc364SMatthew G. Knepley .seealso: DMPlexTransformType, DMRefine(), DMPlexGetTransformType(), DMPlexSetRefinementUniform() 118*012bc364SMatthew G. Knepley @*/ 119*012bc364SMatthew G. Knepley PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type) 120*012bc364SMatthew G. Knepley { 121*012bc364SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 122*012bc364SMatthew G. Knepley PetscErrorCode ierr; 123*012bc364SMatthew G. Knepley 124*012bc364SMatthew G. Knepley PetscFunctionBegin; 125*012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 126*012bc364SMatthew G. Knepley if (type) PetscValidCharPointer(type, 2); 127*012bc364SMatthew G. Knepley ierr = PetscFree(mesh->transformType);CHKERRQ(ierr); 128*012bc364SMatthew G. Knepley ierr = PetscStrallocpy(type, &mesh->transformType);CHKERRQ(ierr); 129*012bc364SMatthew G. Knepley PetscFunctionReturn(0); 130*012bc364SMatthew G. Knepley } 131*012bc364SMatthew G. Knepley 132*012bc364SMatthew G. Knepley /*@C 133*012bc364SMatthew G. Knepley DMPlexGetTransformType - Retrieve the transform type for uniform refinement 134*012bc364SMatthew G. Knepley 135*012bc364SMatthew G. Knepley Input Parameter: 136*012bc364SMatthew G. Knepley . dm - The DM 137*012bc364SMatthew G. Knepley 138*012bc364SMatthew G. Knepley Output Parameter: 139*012bc364SMatthew G. Knepley . type - The transform type for uniform refinement 140*012bc364SMatthew G. Knepley 141*012bc364SMatthew G. Knepley Level: developer 142*012bc364SMatthew G. Knepley 143*012bc364SMatthew G. Knepley .seealso: DMPlexTransformType, DMRefine(), DMPlexSetTransformType(), DMPlexGetRefinementUniform() 144*012bc364SMatthew G. Knepley @*/ 145*012bc364SMatthew G. Knepley PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type) 146*012bc364SMatthew G. Knepley { 147*012bc364SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 148*012bc364SMatthew G. Knepley 149*012bc364SMatthew G. Knepley PetscFunctionBegin; 150*012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 151*012bc364SMatthew G. Knepley PetscValidPointer(type, 2); 152*012bc364SMatthew G. Knepley *type = mesh->transformType; 153*012bc364SMatthew G. Knepley PetscFunctionReturn(0); 154*012bc364SMatthew G. Knepley } 155*012bc364SMatthew G. Knepley 1560e2b6761SMatthew G. Knepley /*@ 1570e2b6761SMatthew G. Knepley DMPlexSetRefinementUniform - Set the flag for uniform refinement 1580e2b6761SMatthew G. Knepley 1590e2b6761SMatthew G. Knepley Input Parameters: 1600e2b6761SMatthew G. Knepley + dm - The DM 1610e2b6761SMatthew G. Knepley - refinementUniform - The flag for uniform refinement 1620e2b6761SMatthew G. Knepley 1630e2b6761SMatthew G. Knepley Level: developer 1640e2b6761SMatthew G. Knepley 1650e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexGetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit() 1660e2b6761SMatthew G. Knepley @*/ 16775d3a19aSMatthew G. Knepley PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform) 16875d3a19aSMatthew G. Knepley { 16975d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 17075d3a19aSMatthew G. Knepley 17175d3a19aSMatthew G. Knepley PetscFunctionBegin; 172*012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 17375d3a19aSMatthew G. Knepley mesh->refinementUniform = refinementUniform; 17475d3a19aSMatthew G. Knepley PetscFunctionReturn(0); 17575d3a19aSMatthew G. Knepley } 17675d3a19aSMatthew G. Knepley 1770e2b6761SMatthew G. Knepley /*@ 1780e2b6761SMatthew G. Knepley DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement 1790e2b6761SMatthew G. Knepley 1800e2b6761SMatthew G. Knepley Input Parameter: 1810e2b6761SMatthew G. Knepley . dm - The DM 1820e2b6761SMatthew G. Knepley 1830e2b6761SMatthew G. Knepley Output Parameter: 1840e2b6761SMatthew G. Knepley . refinementUniform - The flag for uniform refinement 1850e2b6761SMatthew G. Knepley 1860e2b6761SMatthew G. Knepley Level: developer 1870e2b6761SMatthew G. Knepley 1880e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexSetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit() 1890e2b6761SMatthew G. Knepley @*/ 19075d3a19aSMatthew G. Knepley PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform) 19175d3a19aSMatthew G. Knepley { 19275d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 19375d3a19aSMatthew G. Knepley 19475d3a19aSMatthew G. Knepley PetscFunctionBegin; 195*012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 19675d3a19aSMatthew G. Knepley PetscValidPointer(refinementUniform, 2); 19775d3a19aSMatthew G. Knepley *refinementUniform = mesh->refinementUniform; 19875d3a19aSMatthew G. Knepley PetscFunctionReturn(0); 19975d3a19aSMatthew G. Knepley } 20075d3a19aSMatthew G. Knepley 2010e2b6761SMatthew G. Knepley /*@ 2020e2b6761SMatthew G. Knepley DMPlexSetRefinementLimit - Set the maximum cell volume for refinement 2030e2b6761SMatthew G. Knepley 2040e2b6761SMatthew G. Knepley Input Parameters: 2050e2b6761SMatthew G. Knepley + dm - The DM 2060e2b6761SMatthew G. Knepley - refinementLimit - The maximum cell volume in the refined mesh 2070e2b6761SMatthew G. Knepley 2080e2b6761SMatthew G. Knepley Level: developer 2090e2b6761SMatthew G. Knepley 2100e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexGetRefinementLimit(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform() 2110e2b6761SMatthew G. Knepley @*/ 21275d3a19aSMatthew G. Knepley PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit) 21375d3a19aSMatthew G. Knepley { 21475d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 21575d3a19aSMatthew G. Knepley 21675d3a19aSMatthew G. Knepley PetscFunctionBegin; 21775d3a19aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 21875d3a19aSMatthew G. Knepley mesh->refinementLimit = refinementLimit; 21975d3a19aSMatthew G. Knepley PetscFunctionReturn(0); 22075d3a19aSMatthew G. Knepley } 22175d3a19aSMatthew G. Knepley 2220e2b6761SMatthew G. Knepley /*@ 2230e2b6761SMatthew G. Knepley DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement 2240e2b6761SMatthew G. Knepley 2250e2b6761SMatthew G. Knepley Input Parameter: 2260e2b6761SMatthew G. Knepley . dm - The DM 2270e2b6761SMatthew G. Knepley 2280e2b6761SMatthew G. Knepley Output Parameter: 2290e2b6761SMatthew G. Knepley . refinementLimit - The maximum cell volume in the refined mesh 2300e2b6761SMatthew G. Knepley 2310e2b6761SMatthew G. Knepley Level: developer 2320e2b6761SMatthew G. Knepley 2330e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexSetRefinementLimit(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform() 2340e2b6761SMatthew G. Knepley @*/ 23575d3a19aSMatthew G. Knepley PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit) 23675d3a19aSMatthew G. Knepley { 23775d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 23875d3a19aSMatthew G. Knepley 23975d3a19aSMatthew G. Knepley PetscFunctionBegin; 24075d3a19aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 24175d3a19aSMatthew G. Knepley PetscValidPointer(refinementLimit, 2); 24275d3a19aSMatthew G. Knepley /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */ 24375d3a19aSMatthew G. Knepley *refinementLimit = mesh->refinementLimit; 24475d3a19aSMatthew G. Knepley PetscFunctionReturn(0); 24575d3a19aSMatthew G. Knepley } 24675d3a19aSMatthew G. Knepley 247b28003e6SMatthew G. Knepley /*@ 248b28003e6SMatthew G. Knepley DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement 249b28003e6SMatthew G. Knepley 250b28003e6SMatthew G. Knepley Input Parameters: 251b28003e6SMatthew G. Knepley + dm - The DM 252b28003e6SMatthew G. Knepley - refinementFunc - Function giving the maximum cell volume in the refined mesh 253b28003e6SMatthew G. Knepley 254b28003e6SMatthew G. Knepley Note: The calling sequence is refinementFunc(coords, limit) 255b28003e6SMatthew G. Knepley $ coords - Coordinates of the current point, usually a cell centroid 256b28003e6SMatthew G. Knepley $ limit - The maximum cell volume for a cell containing this point 257b28003e6SMatthew G. Knepley 258b28003e6SMatthew G. Knepley Level: developer 259b28003e6SMatthew G. Knepley 260b28003e6SMatthew G. Knepley .seealso: DMRefine(), DMPlexGetRefinementFunction(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit() 261b28003e6SMatthew G. Knepley @*/ 262b28003e6SMatthew G. Knepley PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *)) 263b28003e6SMatthew G. Knepley { 264b28003e6SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 265b28003e6SMatthew G. Knepley 266b28003e6SMatthew G. Knepley PetscFunctionBegin; 267b28003e6SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 268b28003e6SMatthew G. Knepley mesh->refinementFunc = refinementFunc; 269b28003e6SMatthew G. Knepley PetscFunctionReturn(0); 270b28003e6SMatthew G. Knepley } 271b28003e6SMatthew G. Knepley 272b28003e6SMatthew G. Knepley /*@ 273b28003e6SMatthew G. Knepley DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement 274b28003e6SMatthew G. Knepley 275b28003e6SMatthew G. Knepley Input Parameter: 276b28003e6SMatthew G. Knepley . dm - The DM 277b28003e6SMatthew G. Knepley 278b28003e6SMatthew G. Knepley Output Parameter: 279b28003e6SMatthew G. Knepley . refinementFunc - Function giving the maximum cell volume in the refined mesh 280b28003e6SMatthew G. Knepley 281b28003e6SMatthew G. Knepley Note: The calling sequence is refinementFunc(coords, limit) 282b28003e6SMatthew G. Knepley $ coords - Coordinates of the current point, usually a cell centroid 283b28003e6SMatthew G. Knepley $ limit - The maximum cell volume for a cell containing this point 284b28003e6SMatthew G. Knepley 285b28003e6SMatthew G. Knepley Level: developer 286b28003e6SMatthew G. Knepley 287b28003e6SMatthew G. Knepley .seealso: DMRefine(), DMPlexSetRefinementFunction(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit() 288b28003e6SMatthew G. Knepley @*/ 289b28003e6SMatthew G. Knepley PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal [], PetscReal *)) 290b28003e6SMatthew G. Knepley { 291b28003e6SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 292b28003e6SMatthew G. Knepley 293b28003e6SMatthew G. Knepley PetscFunctionBegin; 294b28003e6SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 295b28003e6SMatthew G. Knepley PetscValidPointer(refinementFunc, 2); 296b28003e6SMatthew G. Knepley *refinementFunc = mesh->refinementFunc; 297b28003e6SMatthew G. Knepley PetscFunctionReturn(0); 298b28003e6SMatthew G. Knepley } 299b28003e6SMatthew G. Knepley 300*012bc364SMatthew G. Knepley PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm) 3010d1cd5e0SMatthew G. Knepley { 302492b8470SStefano Zampini PetscBool isUniform; 3030d1cd5e0SMatthew G. Knepley PetscErrorCode ierr; 3040d1cd5e0SMatthew G. Knepley 3050d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 3060d1cd5e0SMatthew G. Knepley ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr); 307412e9a14SMatthew G. Knepley ierr = DMViewFromOptions(dm, NULL, "-initref_dm_view");CHKERRQ(ierr); 3080d1cd5e0SMatthew G. Knepley if (isUniform) { 309*012bc364SMatthew G. Knepley DMPlexTransform tr; 31051a74b61SMatthew G. Knepley DM cdm, rcdm; 311*012bc364SMatthew G. Knepley DMPlexTransformType trType; 312*012bc364SMatthew G. Knepley const char *prefix; 313*012bc364SMatthew G. Knepley PetscOptions options; 3140d1cd5e0SMatthew G. Knepley 315*012bc364SMatthew G. Knepley ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr);CHKERRQ(ierr); 316*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetDM(tr, dm);CHKERRQ(ierr); 317*012bc364SMatthew G. Knepley ierr = DMPlexGetTransformType(dm, &trType);CHKERRQ(ierr); 318*012bc364SMatthew G. Knepley if (trType) {ierr = DMPlexTransformSetType(tr, trType);CHKERRQ(ierr);} 319*012bc364SMatthew G. Knepley ierr = PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix);CHKERRQ(ierr); 320*012bc364SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) tr, prefix);CHKERRQ(ierr); 321*012bc364SMatthew G. Knepley ierr = PetscObjectGetOptions((PetscObject) dm, &options);CHKERRQ(ierr); 322*012bc364SMatthew G. Knepley ierr = PetscObjectSetOptions((PetscObject) tr, options);CHKERRQ(ierr); 323*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetFromOptions(tr);CHKERRQ(ierr); 324*012bc364SMatthew G. Knepley ierr = PetscObjectSetOptions((PetscObject) tr, NULL);CHKERRQ(ierr); 325*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetUp(tr);CHKERRQ(ierr); 326*012bc364SMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) tr, NULL, "-dm_plex_transform_view");CHKERRQ(ierr); 327*012bc364SMatthew G. Knepley ierr = DMPlexTransformApply(tr, dm, rdm);CHKERRQ(ierr); 328*012bc364SMatthew G. Knepley ierr = DMPlexSetRegularRefinement(*rdm, PETSC_TRUE);CHKERRQ(ierr); 329*012bc364SMatthew G. Knepley ierr = DMCopyDisc(dm, *rdm);CHKERRQ(ierr); 33051a74b61SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 331*012bc364SMatthew G. Knepley ierr = DMGetCoordinateDM(*rdm, &rcdm);CHKERRQ(ierr); 33251a74b61SMatthew G. Knepley ierr = DMCopyDisc(cdm, rcdm);CHKERRQ(ierr); 333*012bc364SMatthew G. Knepley ierr = DMPlexTransformCreateDiscLabels(tr, *rdm);CHKERRQ(ierr); 334*012bc364SMatthew G. Knepley ierr = DMPlexTransformDestroy(&tr);CHKERRQ(ierr); 3350d1cd5e0SMatthew G. Knepley } else { 336*012bc364SMatthew G. Knepley ierr = DMPlexRefine_Internal(dm, NULL, rdm);CHKERRQ(ierr); 3370d1cd5e0SMatthew G. Knepley } 338*012bc364SMatthew G. Knepley if (*rdm) { 339*012bc364SMatthew G. Knepley ((DM_Plex *) (*rdm)->data)->printFEM = ((DM_Plex *) dm->data)->printFEM; 340*012bc364SMatthew G. Knepley ((DM_Plex *) (*rdm)->data)->printL2 = ((DM_Plex *) dm->data)->printL2; 3416d7c9049SMatthew G. Knepley } 342*012bc364SMatthew G. Knepley ierr = DMViewFromOptions(dm, NULL, "-postref_dm_view");CHKERRQ(ierr); 3430d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 3440d1cd5e0SMatthew G. Knepley } 3450d1cd5e0SMatthew G. Knepley 346*012bc364SMatthew G. Knepley PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[]) 3470d1cd5e0SMatthew G. Knepley { 3480d1cd5e0SMatthew G. Knepley DM cdm = dm; 3490d1cd5e0SMatthew G. Knepley PetscInt r; 3500d1cd5e0SMatthew G. Knepley PetscBool isUniform, localized; 3510d1cd5e0SMatthew G. Knepley PetscErrorCode ierr; 3520d1cd5e0SMatthew G. Knepley 3530d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 3540d1cd5e0SMatthew G. Knepley ierr = DMPlexGetRefinementUniform(dm, &isUniform);CHKERRQ(ierr); 3550d1cd5e0SMatthew G. Knepley ierr = DMGetCoordinatesLocalized(dm, &localized);CHKERRQ(ierr); 3560d1cd5e0SMatthew G. Knepley if (isUniform) { 3570d1cd5e0SMatthew G. Knepley for (r = 0; r < nlevels; ++r) { 358*012bc364SMatthew G. Knepley DMPlexTransform tr; 35951a74b61SMatthew G. Knepley DM codm, rcodm; 360*012bc364SMatthew G. Knepley const char *prefix; 3610d1cd5e0SMatthew G. Knepley 362*012bc364SMatthew G. Knepley ierr = DMPlexTransformCreate(PetscObjectComm((PetscObject) cdm), &tr);CHKERRQ(ierr); 363*012bc364SMatthew G. Knepley ierr = PetscObjectGetOptionsPrefix((PetscObject) cdm, &prefix);CHKERRQ(ierr); 364*012bc364SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) tr, prefix);CHKERRQ(ierr); 365*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetDM(tr, cdm);CHKERRQ(ierr); 366*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetFromOptions(tr);CHKERRQ(ierr); 367*012bc364SMatthew G. Knepley ierr = DMPlexTransformSetUp(tr);CHKERRQ(ierr); 368*012bc364SMatthew G. Knepley ierr = DMPlexTransformApply(tr, cdm, &rdm[r]);CHKERRQ(ierr); 369*012bc364SMatthew G. Knepley ierr = DMSetCoarsenLevel(rdm[r], cdm->leveldown);CHKERRQ(ierr); 370*012bc364SMatthew G. Knepley ierr = DMSetRefineLevel(rdm[r], cdm->levelup+1);CHKERRQ(ierr); 371*012bc364SMatthew G. Knepley ierr = DMCopyDisc(cdm, rdm[r]);CHKERRQ(ierr); 37251a74b61SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &codm);CHKERRQ(ierr); 373*012bc364SMatthew G. Knepley ierr = DMGetCoordinateDM(rdm[r], &rcodm);CHKERRQ(ierr); 37451a74b61SMatthew G. Knepley ierr = DMCopyDisc(codm, rcodm);CHKERRQ(ierr); 375*012bc364SMatthew G. Knepley ierr = DMPlexTransformCreateDiscLabels(tr, rdm[r]);CHKERRQ(ierr); 376*012bc364SMatthew G. Knepley ierr = DMSetCoarseDM(rdm[r], cdm);CHKERRQ(ierr); 377*012bc364SMatthew G. Knepley ierr = DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE);CHKERRQ(ierr); 378*012bc364SMatthew G. Knepley if (rdm[r]) { 379*012bc364SMatthew G. Knepley ((DM_Plex *) (rdm[r])->data)->printFEM = ((DM_Plex *) dm->data)->printFEM; 380*012bc364SMatthew G. Knepley ((DM_Plex *) (rdm[r])->data)->printL2 = ((DM_Plex *) dm->data)->printL2; 3816d7c9049SMatthew G. Knepley } 382*012bc364SMatthew G. Knepley cdm = rdm[r]; 383*012bc364SMatthew G. Knepley ierr = DMPlexTransformDestroy(&tr);CHKERRQ(ierr); 3840d1cd5e0SMatthew G. Knepley } 3850d1cd5e0SMatthew G. Knepley } else { 3860d1cd5e0SMatthew G. Knepley for (r = 0; r < nlevels; ++r) { 387*012bc364SMatthew G. Knepley ierr = DMRefine(cdm, PetscObjectComm((PetscObject) dm), &rdm[r]);CHKERRQ(ierr); 388*012bc364SMatthew G. Knepley ierr = DMCopyDisc(cdm, rdm[r]);CHKERRQ(ierr); 389*012bc364SMatthew G. Knepley if (localized) {ierr = DMLocalizeCoordinates(rdm[r]);CHKERRQ(ierr);} 390*012bc364SMatthew G. Knepley ierr = DMSetCoarseDM(rdm[r], cdm);CHKERRQ(ierr); 391*012bc364SMatthew G. Knepley if (rdm[r]) { 392*012bc364SMatthew G. Knepley ((DM_Plex *) (rdm[r])->data)->printFEM = ((DM_Plex *) dm->data)->printFEM; 393*012bc364SMatthew G. Knepley ((DM_Plex *) (rdm[r])->data)->printL2 = ((DM_Plex *) dm->data)->printL2; 3946d7c9049SMatthew G. Knepley } 395*012bc364SMatthew G. Knepley cdm = rdm[r]; 3960d1cd5e0SMatthew G. Knepley } 3970d1cd5e0SMatthew G. Knepley } 3980d1cd5e0SMatthew G. Knepley PetscFunctionReturn(0); 3990d1cd5e0SMatthew G. Knepley } 400