1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2412e9a14SMatthew G. Knepley #include <petsc/private/petscfeimpl.h> /* For PetscFEInterpolate_Static() */ 3012bc364SMatthew G. Knepley 4cc4c1da9SBarry Smith #include <petscdmplextransform.h> /*I "petscdmplextransform.h" I*/ 575d3a19aSMatthew G. Knepley #include <petscsf.h> 675d3a19aSMatthew G. Knepley 7963fc26aSMatthew G. Knepley /*@ 8a1cb98faSBarry Smith DMPlexCreateProcessSF - Create an `PetscSF` which just has process connectivity 9963fc26aSMatthew G. Knepley 1020f4b53cSBarry Smith Collective 11963fc26aSMatthew G. Knepley 12963fc26aSMatthew G. Knepley Input Parameters: 13a1cb98faSBarry Smith + dm - The `DM` 14a1cb98faSBarry Smith - sfPoint - The `PetscSF` which encodes point connectivity 15963fc26aSMatthew G. Knepley 16963fc26aSMatthew G. Knepley Output Parameters: 1720f4b53cSBarry Smith + processRanks - A list of process neighbors, or `NULL` 1820f4b53cSBarry Smith - sfProcess - An `PetscSF` encoding the process connectivity, or `NULL` 19963fc26aSMatthew G. Knepley 20963fc26aSMatthew G. Knepley Level: developer 21963fc26aSMatthew G. Knepley 221cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSFCreate()`, `DMPlexCreateTwoSidedProcessSF()` 23963fc26aSMatthew G. Knepley @*/ 24d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess) 25d71ae5a4SJacob Faibussowitsch { 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; 316497c311SBarry Smith PetscMPIInt *ranks; 326497c311SBarry Smith PetscInt *ranksNew; 339852e123SBarry Smith PetscMPIInt size; 3475d3a19aSMatthew G. Knepley 3575d3a19aSMatthew G. Knepley PetscFunctionBegin; 36963fc26aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 37963fc26aSMatthew G. Knepley PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 384f572ea9SToby Isaac if (processRanks) PetscAssertPointer(processRanks, 3); 394f572ea9SToby Isaac if (sfProcess) PetscAssertPointer(sfProcess, 4); 409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 419566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints)); 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeaves, &ranks)); 436497c311SBarry Smith for (l = 0; l < numLeaves; ++l) ranks[l] = (PetscMPIInt)remotePoints[l].rank; 446497c311SBarry Smith PetscCall(PetscSortRemoveDupsMPIInt(&numLeaves, ranks)); 459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeaves, &ranksNew)); 469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeaves, &localPointsNew)); 479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numLeaves, &remotePointsNew)); 4875d3a19aSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) { 4975d3a19aSMatthew G. Knepley ranksNew[l] = ranks[l]; 5075d3a19aSMatthew G. Knepley localPointsNew[l] = l; 5175d3a19aSMatthew G. Knepley remotePointsNew[l].index = 0; 526497c311SBarry Smith remotePointsNew[l].rank = (PetscMPIInt)ranksNew[l]; 5375d3a19aSMatthew G. Knepley } 549566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks)); 559566063dSJacob Faibussowitsch if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks)); 569566063dSJacob Faibussowitsch else PetscCall(PetscFree(ranksNew)); 57963fc26aSMatthew G. Knepley if (sfProcess) { 589566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Process SF")); 609566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(*sfProcess)); 619566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER)); 62963fc26aSMatthew G. Knepley } 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6475d3a19aSMatthew G. Knepley } 6575d3a19aSMatthew G. Knepley 662389894bSMatthew G. Knepley /*@ 67a1cb98faSBarry Smith DMPlexCreateCoarsePointIS - Creates an `IS` covering the coarse `DM` chart with the fine points as data 682389894bSMatthew G. Knepley 6920f4b53cSBarry Smith Collective 70012bc364SMatthew G. Knepley 712389894bSMatthew G. Knepley Input Parameter: 72a1cb98faSBarry Smith . dm - The coarse `DM` 732389894bSMatthew G. Knepley 742389894bSMatthew G. Knepley Output Parameter: 75a1cb98faSBarry Smith . fpointIS - The `IS` of all the fine points which exist in the original coarse mesh 762389894bSMatthew G. Knepley 772389894bSMatthew G. Knepley Level: developer 782389894bSMatthew G. Knepley 791cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `IS`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetSubpointIS()` 802389894bSMatthew G. Knepley @*/ 81d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS) 82d71ae5a4SJacob Faibussowitsch { 83012bc364SMatthew G. Knepley DMPlexTransform tr; 84412e9a14SMatthew G. Knepley PetscInt *fpoints; 85327c2912SStefano Zampini PetscInt pStart, pEnd, p, vStart, vEnd, v; 862389894bSMatthew G. Knepley 872389894bSMatthew G. Knepley PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 899566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 909566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 919566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(tr)); 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pEnd - pStart, &fpoints)); 932389894bSMatthew G. Knepley for (p = 0; p < pEnd - pStart; ++p) fpoints[p] = -1; 94412e9a14SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 95012bc364SMatthew G. Knepley PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */ 96327c2912SStefano Zampini 979566063dSJacob Faibussowitsch PetscCall(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew)); 98412e9a14SMatthew G. Knepley fpoints[v - pStart] = vNew; 992389894bSMatthew G. Knepley } 1009566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 1019566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, fpoints, PETSC_OWN_POINTER, fpointIS)); 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1032389894bSMatthew G. Knepley } 1042389894bSMatthew G. Knepley 105012bc364SMatthew G. Knepley /*@C 106012bc364SMatthew G. Knepley DMPlexSetTransformType - Set the transform type for uniform refinement 107012bc364SMatthew G. Knepley 108012bc364SMatthew G. Knepley Input Parameters: 109a1cb98faSBarry Smith + dm - The `DM` 110012bc364SMatthew G. Knepley - type - The transform type for uniform refinement 111012bc364SMatthew G. Knepley 112012bc364SMatthew G. Knepley Level: developer 113012bc364SMatthew G. Knepley 1141cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexGetTransformType()`, `DMPlexSetRefinementUniform()` 115012bc364SMatthew G. Knepley @*/ 116d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type) 117d71ae5a4SJacob Faibussowitsch { 118012bc364SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 119012bc364SMatthew G. Knepley 120012bc364SMatthew G. Knepley PetscFunctionBegin; 121012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 1224f572ea9SToby Isaac if (type) PetscAssertPointer(type, 2); 1239566063dSJacob Faibussowitsch PetscCall(PetscFree(mesh->transformType)); 1249566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(type, &mesh->transformType)); 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126012bc364SMatthew G. Knepley } 127012bc364SMatthew G. Knepley 128012bc364SMatthew G. Knepley /*@C 129012bc364SMatthew G. Knepley DMPlexGetTransformType - Retrieve the transform type for uniform refinement 130012bc364SMatthew G. Knepley 131012bc364SMatthew G. Knepley Input Parameter: 132a1cb98faSBarry Smith . dm - The `DM` 133012bc364SMatthew G. Knepley 134012bc364SMatthew G. Knepley Output Parameter: 135012bc364SMatthew G. Knepley . type - The transform type for uniform refinement 136012bc364SMatthew G. Knepley 137012bc364SMatthew G. Knepley Level: developer 138012bc364SMatthew G. Knepley 1391cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexSetTransformType()`, `DMPlexGetRefinementUniform()` 140012bc364SMatthew G. Knepley @*/ 141d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type) 142d71ae5a4SJacob Faibussowitsch { 143012bc364SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 144012bc364SMatthew G. Knepley 145012bc364SMatthew G. Knepley PetscFunctionBegin; 146012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 1474f572ea9SToby Isaac PetscAssertPointer(type, 2); 148012bc364SMatthew G. Knepley *type = mesh->transformType; 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150012bc364SMatthew G. Knepley } 151012bc364SMatthew G. Knepley 1520e2b6761SMatthew G. Knepley /*@ 1530e2b6761SMatthew G. Knepley DMPlexSetRefinementUniform - Set the flag for uniform refinement 1540e2b6761SMatthew G. Knepley 1550e2b6761SMatthew G. Knepley Input Parameters: 156a1cb98faSBarry Smith + dm - The `DM` 1570e2b6761SMatthew G. Knepley - refinementUniform - The flag for uniform refinement 1580e2b6761SMatthew G. Knepley 1590e2b6761SMatthew G. Knepley Level: developer 1600e2b6761SMatthew G. Knepley 1611cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 1620e2b6761SMatthew G. Knepley @*/ 163d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform) 164d71ae5a4SJacob Faibussowitsch { 16575d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 16675d3a19aSMatthew G. Knepley 16775d3a19aSMatthew G. Knepley PetscFunctionBegin; 168012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 16975d3a19aSMatthew G. Knepley mesh->refinementUniform = refinementUniform; 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17175d3a19aSMatthew G. Knepley } 17275d3a19aSMatthew G. Knepley 1730e2b6761SMatthew G. Knepley /*@ 1740e2b6761SMatthew G. Knepley DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement 1750e2b6761SMatthew G. Knepley 1760e2b6761SMatthew G. Knepley Input Parameter: 177a1cb98faSBarry Smith . dm - The `DM` 1780e2b6761SMatthew G. Knepley 1790e2b6761SMatthew G. Knepley Output Parameter: 1800e2b6761SMatthew G. Knepley . refinementUniform - The flag for uniform refinement 1810e2b6761SMatthew G. Knepley 1820e2b6761SMatthew G. Knepley Level: developer 1830e2b6761SMatthew G. Knepley 1841cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 1850e2b6761SMatthew G. Knepley @*/ 186d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform) 187d71ae5a4SJacob Faibussowitsch { 18875d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 18975d3a19aSMatthew G. Knepley 19075d3a19aSMatthew G. Knepley PetscFunctionBegin; 191012bc364SMatthew G. Knepley PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX); 1924f572ea9SToby Isaac PetscAssertPointer(refinementUniform, 2); 19375d3a19aSMatthew G. Knepley *refinementUniform = mesh->refinementUniform; 1943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19575d3a19aSMatthew G. Knepley } 19675d3a19aSMatthew G. Knepley 1970e2b6761SMatthew G. Knepley /*@ 1980e2b6761SMatthew G. Knepley DMPlexSetRefinementLimit - Set the maximum cell volume for refinement 1990e2b6761SMatthew G. Knepley 2000e2b6761SMatthew G. Knepley Input Parameters: 201a1cb98faSBarry Smith + dm - The `DM` 2020e2b6761SMatthew G. Knepley - refinementLimit - The maximum cell volume in the refined mesh 2030e2b6761SMatthew G. Knepley 2040e2b6761SMatthew G. Knepley Level: developer 2050e2b6761SMatthew G. Knepley 2061cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()` 2070e2b6761SMatthew G. Knepley @*/ 208d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit) 209d71ae5a4SJacob Faibussowitsch { 21075d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 21175d3a19aSMatthew G. Knepley 21275d3a19aSMatthew G. Knepley PetscFunctionBegin; 21375d3a19aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 21475d3a19aSMatthew G. Knepley mesh->refinementLimit = refinementLimit; 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21675d3a19aSMatthew G. Knepley } 21775d3a19aSMatthew G. Knepley 2180e2b6761SMatthew G. Knepley /*@ 2190e2b6761SMatthew G. Knepley DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement 2200e2b6761SMatthew G. Knepley 2210e2b6761SMatthew G. Knepley Input Parameter: 222a1cb98faSBarry Smith . dm - The `DM` 2230e2b6761SMatthew G. Knepley 2240e2b6761SMatthew G. Knepley Output Parameter: 2250e2b6761SMatthew G. Knepley . refinementLimit - The maximum cell volume in the refined mesh 2260e2b6761SMatthew G. Knepley 2270e2b6761SMatthew G. Knepley Level: developer 2280e2b6761SMatthew G. Knepley 2291cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()` 2300e2b6761SMatthew G. Knepley @*/ 231d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit) 232d71ae5a4SJacob Faibussowitsch { 23375d3a19aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 23475d3a19aSMatthew G. Knepley 23575d3a19aSMatthew G. Knepley PetscFunctionBegin; 23675d3a19aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2374f572ea9SToby Isaac PetscAssertPointer(refinementLimit, 2); 23875d3a19aSMatthew G. Knepley /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */ 23975d3a19aSMatthew G. Knepley *refinementLimit = mesh->refinementLimit; 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24175d3a19aSMatthew G. Knepley } 24275d3a19aSMatthew G. Knepley 24360225df5SJacob Faibussowitsch /*@C 244b28003e6SMatthew G. Knepley DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement 245b28003e6SMatthew G. Knepley 246b28003e6SMatthew G. Knepley Input Parameters: 247a1cb98faSBarry Smith + dm - The `DM` 248b28003e6SMatthew G. Knepley - refinementFunc - Function giving the maximum cell volume in the refined mesh 249b28003e6SMatthew G. Knepley 25020f4b53cSBarry Smith Calling Sequence of `refinementFunc`: 251a1cb98faSBarry Smith + coords - Coordinates of the current point, usually a cell centroid 252a1cb98faSBarry Smith - limit - The maximum cell volume for a cell containing this point 253b28003e6SMatthew G. Knepley 254b28003e6SMatthew G. Knepley Level: developer 255b28003e6SMatthew G. Knepley 2561cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 257b28003e6SMatthew G. Knepley @*/ 258a4e35b19SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal coords[], PetscReal *limit)) 259d71ae5a4SJacob Faibussowitsch { 260b28003e6SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 261b28003e6SMatthew G. Knepley 262b28003e6SMatthew G. Knepley PetscFunctionBegin; 263b28003e6SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 264b28003e6SMatthew G. Knepley mesh->refinementFunc = refinementFunc; 2653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 266b28003e6SMatthew G. Knepley } 267b28003e6SMatthew G. Knepley 26860225df5SJacob Faibussowitsch /*@C 269b28003e6SMatthew G. Knepley DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement 270b28003e6SMatthew G. Knepley 271b28003e6SMatthew G. Knepley Input Parameter: 27220f4b53cSBarry Smith . dm - The `DM` 273b28003e6SMatthew G. Knepley 274b28003e6SMatthew G. Knepley Output Parameter: 275b28003e6SMatthew G. Knepley . refinementFunc - Function giving the maximum cell volume in the refined mesh 276b28003e6SMatthew G. Knepley 27720f4b53cSBarry Smith Calling Sequence of `refinementFunc`: 278a1cb98faSBarry Smith + coords - Coordinates of the current point, usually a cell centroid 279a1cb98faSBarry Smith - limit - The maximum cell volume for a cell containing this point 280b28003e6SMatthew G. Knepley 281b28003e6SMatthew G. Knepley Level: developer 282b28003e6SMatthew G. Knepley 2831cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()` 284b28003e6SMatthew G. Knepley @*/ 285a4e35b19SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal coords[], PetscReal *limit)) 286d71ae5a4SJacob Faibussowitsch { 287b28003e6SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 288b28003e6SMatthew G. Knepley 289b28003e6SMatthew G. Knepley PetscFunctionBegin; 290b28003e6SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2914f572ea9SToby Isaac PetscAssertPointer(refinementFunc, 2); 292b28003e6SMatthew G. Knepley *refinementFunc = mesh->refinementFunc; 2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 294b28003e6SMatthew G. Knepley } 295b28003e6SMatthew G. Knepley 296d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm) 297d71ae5a4SJacob Faibussowitsch { 298492b8470SStefano Zampini PetscBool isUniform; 2990d1cd5e0SMatthew G. Knepley 3000d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 3019566063dSJacob Faibussowitsch PetscCall(DMPlexGetRefinementUniform(dm, &isUniform)); 3029566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-initref_dm_view")); 3030d1cd5e0SMatthew G. Knepley if (isUniform) { 304012bc364SMatthew G. Knepley DMPlexTransform tr; 30551a74b61SMatthew G. Knepley DM cdm, rcdm; 306012bc364SMatthew G. Knepley DMPlexTransformType trType; 307012bc364SMatthew G. Knepley const char *prefix; 308012bc364SMatthew G. Knepley PetscOptions options; 309e44f6aebSMatthew G. Knepley PetscInt cDegree; 310e01335dfSMatthew G. Knepley PetscBool useCeed, flg; 311e01335dfSMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 3120d1cd5e0SMatthew G. Knepley 3139566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr)); 3149566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetDM(tr, dm)); 3159566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransformType(dm, &trType)); 3169566063dSJacob Faibussowitsch if (trType) PetscCall(DMPlexTransformSetType(tr, trType)); 3179566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 3189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 3199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 3209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptions((PetscObject)tr, options)); 3219566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetFromOptions(tr)); 3229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL)); 323e01335dfSMatthew G. Knepley PetscCall(PetscOptionsGetString(options, prefix, "-dm_plex_transform_active", name, PETSC_MAX_PATH_LEN, &flg)); 324e01335dfSMatthew G. Knepley if (flg) { 325e01335dfSMatthew G. Knepley PetscCall(DMHasLabel(dm, name, &flg)); 326e01335dfSMatthew G. Knepley if (flg) { 327e01335dfSMatthew G. Knepley DMLabel active; 328e01335dfSMatthew G. Knepley 329e01335dfSMatthew G. Knepley PetscCall(DMGetLabel(dm, name, &active)); 330e01335dfSMatthew G. Knepley PetscCall(DMPlexTransformSetActive(tr, active)); 331e01335dfSMatthew G. Knepley } 332e01335dfSMatthew G. Knepley } 3339566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(tr)); 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view")); 3359566063dSJacob Faibussowitsch PetscCall(DMPlexTransformApply(tr, dm, rdm)); 3369566063dSJacob Faibussowitsch PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE)); 337d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 338d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseCeed(*rdm, useCeed)); 339f4f49eeaSPierre Jolivet PetscCall(DMSetMatType(*rdm, dm->mattype)); 3409566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *rdm)); 3419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 3429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*rdm, &rcdm)); 343e44f6aebSMatthew G. Knepley PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree)); 344*38cc584fSMatthew G. Knepley { 345*38cc584fSMatthew G. Knepley PetscDS cds, rcds; 346*38cc584fSMatthew G. Knepley 3474b87bf80SMatthew G. Knepley PetscCall(DMPlexCreateCoordinateSpace(*rdm, cDegree, PETSC_TRUE, NULL)); 348e44f6aebSMatthew G. Knepley PetscCall(DMGetCoordinateDM(*rdm, &rcdm)); 349*38cc584fSMatthew G. Knepley PetscCall(DMGetDS(cdm, &cds)); 350*38cc584fSMatthew G. Knepley PetscCall(DMGetDS(rcdm, &rcds)); 351*38cc584fSMatthew G. Knepley PetscCall(PetscDSCopyConstants(cds, rcds)); 352e44f6aebSMatthew G. Knepley } 353d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(cdm, &useCeed)); 354d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseCeed(rcdm, useCeed)); 355d2b2dc1eSMatthew G. Knepley if (useCeed) { 356d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseMatClosurePermutation(rcdm, PETSC_FALSE)); 357d2b2dc1eSMatthew G. Knepley PetscCall(DMUseTensorOrder(rcdm, PETSC_TRUE)); 358d2b2dc1eSMatthew G. Knepley } 3599566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm)); 3609566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 3610d1cd5e0SMatthew G. Knepley } else { 3629566063dSJacob Faibussowitsch PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm)); 3630d1cd5e0SMatthew G. Knepley } 364012bc364SMatthew G. Knepley if (*rdm) { 365012bc364SMatthew G. Knepley ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 366012bc364SMatthew G. Knepley ((DM_Plex *)(*rdm)->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 3676d7c9049SMatthew G. Knepley } 368e01335dfSMatthew G. Knepley PetscCall(DMViewFromOptions(*rdm, NULL, "-postref_dm_view")); 3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3700d1cd5e0SMatthew G. Knepley } 3710d1cd5e0SMatthew G. Knepley 372d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[]) 373d71ae5a4SJacob Faibussowitsch { 3740d1cd5e0SMatthew G. Knepley DM cdm = dm; 3750d1cd5e0SMatthew G. Knepley PetscInt r; 376d2b2dc1eSMatthew G. Knepley PetscBool isUniform, localized, useCeed; 3770d1cd5e0SMatthew G. Knepley 3780d1cd5e0SMatthew G. Knepley PetscFunctionBegin; 3799566063dSJacob Faibussowitsch PetscCall(DMPlexGetRefinementUniform(dm, &isUniform)); 3809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalized(dm, &localized)); 3810d1cd5e0SMatthew G. Knepley if (isUniform) { 3820d1cd5e0SMatthew G. Knepley for (r = 0; r < nlevels; ++r) { 383012bc364SMatthew G. Knepley DMPlexTransform tr; 38451a74b61SMatthew G. Knepley DM codm, rcodm; 385012bc364SMatthew G. Knepley const char *prefix; 3860d1cd5e0SMatthew G. Knepley 3879566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr)); 3889566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix)); 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix)); 3909566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetDM(tr, cdm)); 3919566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetFromOptions(tr)); 3929566063dSJacob Faibussowitsch PetscCall(DMPlexTransformSetUp(tr)); 3939566063dSJacob Faibussowitsch PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r])); 3949566063dSJacob Faibussowitsch PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown)); 3959566063dSJacob Faibussowitsch PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1)); 396bb393918Smarkadams4 PetscCall(DMSetMatType(rdm[r], dm->mattype)); 397d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 398d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseCeed(rdm[r], useCeed)); 3999566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(cdm, rdm[r])); 4009566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &codm)); 4019566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(rdm[r], &rcodm)); 4029566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(codm, rcodm)); 403d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(codm, &useCeed)); 404d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseCeed(rcodm, useCeed)); 405d2b2dc1eSMatthew G. Knepley if (useCeed) { 406d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseMatClosurePermutation(rcodm, PETSC_FALSE)); 407d2b2dc1eSMatthew G. Knepley PetscCall(DMUseTensorOrder(rcodm, PETSC_TRUE)); 408d2b2dc1eSMatthew G. Knepley } 4099566063dSJacob Faibussowitsch PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r])); 4109566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(rdm[r], cdm)); 4119566063dSJacob Faibussowitsch PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE)); 412012bc364SMatthew G. Knepley if (rdm[r]) { 413f4f49eeaSPierre Jolivet ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 414f4f49eeaSPierre Jolivet ((DM_Plex *)rdm[r]->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 4156d7c9049SMatthew G. Knepley } 416012bc364SMatthew G. Knepley cdm = rdm[r]; 4179566063dSJacob Faibussowitsch PetscCall(DMPlexTransformDestroy(&tr)); 4180d1cd5e0SMatthew G. Knepley } 4190d1cd5e0SMatthew G. Knepley } else { 4200d1cd5e0SMatthew G. Knepley for (r = 0; r < nlevels; ++r) { 4219566063dSJacob Faibussowitsch PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r])); 422d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 423d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexSetUseCeed(rdm[r], useCeed)); 4249566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(cdm, rdm[r])); 4259566063dSJacob Faibussowitsch if (localized) PetscCall(DMLocalizeCoordinates(rdm[r])); 4269566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(rdm[r], cdm)); 427012bc364SMatthew G. Knepley if (rdm[r]) { 428f4f49eeaSPierre Jolivet ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM; 429f4f49eeaSPierre Jolivet ((DM_Plex *)rdm[r]->data)->printL2 = ((DM_Plex *)dm->data)->printL2; 4306d7c9049SMatthew G. Knepley } 431012bc364SMatthew G. Knepley cdm = rdm[r]; 4320d1cd5e0SMatthew G. Knepley } 4330d1cd5e0SMatthew G. Knepley } 4343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4350d1cd5e0SMatthew G. Knepley } 436