16858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 26858538eSMatthew G. Knepley 36858538eSMatthew G. Knepley #include <petscdmplex.h> 46858538eSMatthew G. Knepley 56858538eSMatthew G. Knepley /*@C 66858538eSMatthew G. Knepley DMGetPeriodicity - Get the description of mesh periodicity 76858538eSMatthew G. Knepley 842b02d39SMatthew G. Knepley Not collective 942b02d39SMatthew G. Knepley 106858538eSMatthew G. Knepley Input Parameter: 1120f4b53cSBarry Smith . dm - The `DM` object 126858538eSMatthew G. Knepley 136858538eSMatthew G. Knepley Output Parameters: 146858538eSMatthew G. Knepley + maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates 1520f4b53cSBarry Smith . Lstart - If we assume the mesh is a torus, this is the start of each coordinate, or `NULL` for 0.0 166858538eSMatthew G. Knepley - L - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0 176858538eSMatthew G. Knepley 186858538eSMatthew G. Knepley Level: developer 196858538eSMatthew G. Knepley 2042747ad1SJacob Faibussowitsch .seealso: `DM` 216858538eSMatthew G. Knepley @*/ 225d83a8b1SBarry Smith PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal *maxCell[], const PetscReal *Lstart[], const PetscReal *L[]) 23d71ae5a4SJacob Faibussowitsch { 246858538eSMatthew G. Knepley PetscFunctionBegin; 256858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 266858538eSMatthew G. Knepley if (maxCell) *maxCell = dm->maxCell; 274fb89dddSMatthew G. Knepley if (Lstart) *Lstart = dm->Lstart; 284fb89dddSMatthew G. Knepley if (L) *L = dm->L; 293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 306858538eSMatthew G. Knepley } 316858538eSMatthew G. Knepley 32cc4c1da9SBarry Smith /*@ 336858538eSMatthew G. Knepley DMSetPeriodicity - Set the description of mesh periodicity 346858538eSMatthew G. Knepley 3542b02d39SMatthew G. Knepley Logically Collective 3642b02d39SMatthew G. Knepley 376858538eSMatthew G. Knepley Input Parameters: 3820f4b53cSBarry Smith + dm - The `DM` object 3920f4b53cSBarry Smith . maxCell - Over distances greater than this, we can assume a point has crossed over to another sheet, when trying to localize cell coordinates. Pass `NULL` to remove such information. 4020f4b53cSBarry Smith . Lstart - If we assume the mesh is a torus, this is the start of each coordinate, or `NULL` for 0.0 416858538eSMatthew G. Knepley - L - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0 426858538eSMatthew G. Knepley 436858538eSMatthew G. Knepley Level: developer 446858538eSMatthew G. Knepley 4520f4b53cSBarry Smith .seealso: `DM`, `DMGetPeriodicity()` 466858538eSMatthew G. Knepley @*/ 47d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal Lstart[], const PetscReal L[]) 48d71ae5a4SJacob Faibussowitsch { 496858538eSMatthew G. Knepley PetscInt dim, d; 506858538eSMatthew G. Knepley 516858538eSMatthew G. Knepley PetscFunctionBegin; 526858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 534f572ea9SToby Isaac if (maxCell) PetscAssertPointer(maxCell, 2); 544f572ea9SToby Isaac if (Lstart) PetscAssertPointer(Lstart, 3); 554f572ea9SToby Isaac if (L) PetscAssertPointer(L, 4); 566858538eSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 576858538eSMatthew G. Knepley if (maxCell) { 586858538eSMatthew G. Knepley if (!dm->maxCell) PetscCall(PetscMalloc1(dim, &dm->maxCell)); 596858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) dm->maxCell[d] = maxCell[d]; 606858538eSMatthew G. Knepley } else { /* remove maxCell information to disable automatic computation of localized vertices */ 616858538eSMatthew G. Knepley PetscCall(PetscFree(dm->maxCell)); 626858538eSMatthew G. Knepley dm->maxCell = NULL; 636858538eSMatthew G. Knepley } 644fb89dddSMatthew G. Knepley if (Lstart) { 654fb89dddSMatthew G. Knepley if (!dm->Lstart) PetscCall(PetscMalloc1(dim, &dm->Lstart)); 664fb89dddSMatthew G. Knepley for (d = 0; d < dim; ++d) dm->Lstart[d] = Lstart[d]; 674fb89dddSMatthew G. Knepley } else { /* remove L information to disable automatic computation of localized vertices */ 684fb89dddSMatthew G. Knepley PetscCall(PetscFree(dm->Lstart)); 694fb89dddSMatthew G. Knepley dm->Lstart = NULL; 704fb89dddSMatthew G. Knepley } 716858538eSMatthew G. Knepley if (L) { 726858538eSMatthew G. Knepley if (!dm->L) PetscCall(PetscMalloc1(dim, &dm->L)); 736858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) dm->L[d] = L[d]; 746858538eSMatthew G. Knepley } else { /* remove L information to disable automatic computation of localized vertices */ 756858538eSMatthew G. Knepley PetscCall(PetscFree(dm->L)); 766858538eSMatthew G. Knepley dm->L = NULL; 776858538eSMatthew G. Knepley } 786858538eSMatthew G. Knepley PetscCheck((dm->maxCell && dm->L) || (!dm->maxCell && !dm->L), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot set only one of maxCell/L"); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 806858538eSMatthew G. Knepley } 816858538eSMatthew G. Knepley 826858538eSMatthew G. Knepley /*@ 836858538eSMatthew G. Knepley DMLocalizeCoordinate - If a mesh is periodic (a torus with lengths L_i, some of which can be infinite), project the coordinate onto [0, L_i) in each dimension. 846858538eSMatthew G. Knepley 856858538eSMatthew G. Knepley Input Parameters: 8620f4b53cSBarry Smith + dm - The `DM` 876858538eSMatthew G. Knepley . in - The input coordinate point (dim numbers) 886858538eSMatthew G. Knepley - endpoint - Include the endpoint L_i 896858538eSMatthew G. Knepley 906858538eSMatthew G. Knepley Output Parameter: 910b3275a6SBarry Smith . out - The localized coordinate point (dim numbers) 926858538eSMatthew G. Knepley 936858538eSMatthew G. Knepley Level: developer 946858538eSMatthew G. Knepley 9520f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()` 966858538eSMatthew G. Knepley @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinate(DM dm, const PetscScalar in[], PetscBool endpoint, PetscScalar out[]) 98d71ae5a4SJacob Faibussowitsch { 996858538eSMatthew G. Knepley PetscInt dim, d; 1006858538eSMatthew G. Knepley 1016858538eSMatthew G. Knepley PetscFunctionBegin; 1026858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &dim)); 1036858538eSMatthew G. Knepley if (!dm->maxCell) { 1046858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] = in[d]; 1056858538eSMatthew G. Knepley } else { 1066858538eSMatthew G. Knepley if (endpoint) { 1076858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 1086858538eSMatthew G. Knepley if ((PetscAbsReal(PetscRealPart(in[d]) / dm->L[d] - PetscFloorReal(PetscRealPart(in[d]) / dm->L[d])) < PETSC_SMALL) && (PetscRealPart(in[d]) / dm->L[d] > PETSC_SMALL)) { 1096858538eSMatthew G. Knepley out[d] = in[d] - dm->L[d] * (PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]) - 1); 1106858538eSMatthew G. Knepley } else { 1116858538eSMatthew G. Knepley out[d] = in[d] - dm->L[d] * PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]); 1126858538eSMatthew G. Knepley } 1136858538eSMatthew G. Knepley } 1146858538eSMatthew G. Knepley } else { 115ad540459SPierre Jolivet for (d = 0; d < dim; ++d) out[d] = in[d] - dm->L[d] * PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]); 1166858538eSMatthew G. Knepley } 1176858538eSMatthew G. Knepley } 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1196858538eSMatthew G. Knepley } 1206858538eSMatthew G. Knepley 1216858538eSMatthew G. Knepley /* 1226858538eSMatthew G. Knepley DMLocalizeCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer. 1236858538eSMatthew G. Knepley 1246858538eSMatthew G. Knepley Input Parameters: 12520f4b53cSBarry Smith + dm - The `DM` 1266858538eSMatthew G. Knepley . dim - The spatial dimension 1276858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it 1286858538eSMatthew G. Knepley - in - The input coordinate point (dim numbers) 1296858538eSMatthew G. Knepley 1306858538eSMatthew G. Knepley Output Parameter: 1310b3275a6SBarry Smith . out - The localized coordinate point (dim numbers) 1326858538eSMatthew G. Knepley 1336858538eSMatthew G. Knepley Level: developer 1346858538eSMatthew G. Knepley 13520f4b53cSBarry Smith Note: 13620f4b53cSBarry Smith This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually the one of the vertices on a containing cell 1376858538eSMatthew G. Knepley 13820f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()` 1396858538eSMatthew G. Knepley */ 140d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[]) 141d71ae5a4SJacob Faibussowitsch { 1426858538eSMatthew G. Knepley PetscInt d; 1436858538eSMatthew G. Knepley 1446858538eSMatthew G. Knepley PetscFunctionBegin; 1456858538eSMatthew G. Knepley if (!dm->maxCell) { 1466858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] = in[d]; 1476858538eSMatthew G. Knepley } else { 1486858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 1496858538eSMatthew G. Knepley if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) { 1506858538eSMatthew G. Knepley out[d] = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d]; 1516858538eSMatthew G. Knepley } else { 1526858538eSMatthew G. Knepley out[d] = in[d]; 1536858538eSMatthew G. Knepley } 1546858538eSMatthew G. Knepley } 1556858538eSMatthew G. Knepley } 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1576858538eSMatthew G. Knepley } 1586858538eSMatthew G. Knepley 159d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim, const PetscReal anchor[], const PetscReal in[], PetscReal out[]) 160d71ae5a4SJacob Faibussowitsch { 1616858538eSMatthew G. Knepley PetscInt d; 1626858538eSMatthew G. Knepley 1636858538eSMatthew G. Knepley PetscFunctionBegin; 1646858538eSMatthew G. Knepley if (!dm->maxCell) { 1656858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] = in[d]; 1666858538eSMatthew G. Knepley } else { 1676858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 1686858538eSMatthew G. Knepley if ((dm->L[d] > 0.0) && (PetscAbsReal(anchor[d] - in[d]) > dm->maxCell[d])) { 1696858538eSMatthew G. Knepley out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d]; 1706858538eSMatthew G. Knepley } else { 1716858538eSMatthew G. Knepley out[d] = in[d]; 1726858538eSMatthew G. Knepley } 1736858538eSMatthew G. Knepley } 1746858538eSMatthew G. Knepley } 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1766858538eSMatthew G. Knepley } 1776858538eSMatthew G. Knepley 1786858538eSMatthew G. Knepley /* 1796858538eSMatthew G. Knepley DMLocalizeAddCoordinate_Internal - If a mesh is periodic, and the input point is far from the anchor, pick the coordinate sheet of the torus which moves it closer. 1806858538eSMatthew G. Knepley 1816858538eSMatthew G. Knepley Input Parameters: 18220f4b53cSBarry Smith + dm - The `DM` 1836858538eSMatthew G. Knepley . dim - The spatial dimension 1846858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it 1856858538eSMatthew G. Knepley . in - The input coordinate delta (dim numbers) 1866858538eSMatthew G. Knepley - out - The input coordinate point (dim numbers) 1876858538eSMatthew G. Knepley 1886858538eSMatthew G. Knepley Output Parameter: 1896858538eSMatthew G. Knepley . out - The localized coordinate in + out 1906858538eSMatthew G. Knepley 1916858538eSMatthew G. Knepley Level: developer 1926858538eSMatthew G. Knepley 19320f4b53cSBarry Smith Note: 1940b3275a6SBarry Smith This is meant to get a set of coordinates close to each other, as in a cell. The anchor is usually one of the vertices on a containing cell 1956858538eSMatthew G. Knepley 19620f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMLocalizeCoordinate()` 1976858538eSMatthew G. Knepley */ 198d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeAddCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[]) 199d71ae5a4SJacob Faibussowitsch { 2006858538eSMatthew G. Knepley PetscInt d; 2016858538eSMatthew G. Knepley 2026858538eSMatthew G. Knepley PetscFunctionBegin; 2036858538eSMatthew G. Knepley if (!dm->maxCell) { 2046858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) out[d] += in[d]; 2056858538eSMatthew G. Knepley } else { 2066858538eSMatthew G. Knepley for (d = 0; d < dim; ++d) { 2076858538eSMatthew G. Knepley const PetscReal maxC = dm->maxCell[d]; 2086858538eSMatthew G. Knepley 2096858538eSMatthew G. Knepley if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > maxC)) { 2106858538eSMatthew G. Knepley const PetscScalar newCoord = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d]; 2116858538eSMatthew G. Knepley 2126858538eSMatthew G. Knepley if (PetscAbsScalar(newCoord - anchor[d]) > maxC) 2136858538eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "%" PetscInt_FMT "-Coordinate %g more than %g away from anchor %g", d, (double)PetscRealPart(in[d]), (double)maxC, (double)PetscRealPart(anchor[d])); 2146858538eSMatthew G. Knepley out[d] += newCoord; 2156858538eSMatthew G. Knepley } else { 2166858538eSMatthew G. Knepley out[d] += in[d]; 2176858538eSMatthew G. Knepley } 2186858538eSMatthew G. Knepley } 2196858538eSMatthew G. Knepley } 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2216858538eSMatthew G. Knepley } 2226858538eSMatthew G. Knepley 2236858538eSMatthew G. Knepley /*@ 22420f4b53cSBarry Smith DMGetCoordinatesLocalizedLocal - Check if the `DM` coordinates have been localized for cells on this process 2256858538eSMatthew G. Knepley 22620f4b53cSBarry Smith Not Collective 2276858538eSMatthew G. Knepley 2286858538eSMatthew G. Knepley Input Parameter: 22920f4b53cSBarry Smith . dm - The `DM` 2306858538eSMatthew G. Knepley 2316858538eSMatthew G. Knepley Output Parameter: 232a4e35b19SJacob Faibussowitsch . areLocalized - `PETSC_TRUE` if localized 2336858538eSMatthew G. Knepley 2346858538eSMatthew G. Knepley Level: developer 2356858538eSMatthew G. Knepley 23620f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMGetCoordinatesLocalized()`, `DMSetPeriodicity()` 2376858538eSMatthew G. Knepley @*/ 238d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalizedLocal(DM dm, PetscBool *areLocalized) 239d71ae5a4SJacob Faibussowitsch { 2406858538eSMatthew G. Knepley PetscFunctionBegin; 2416858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2424f572ea9SToby Isaac PetscAssertPointer(areLocalized, 2); 2436858538eSMatthew G. Knepley *areLocalized = dm->coordinates[1].dim < 0 ? PETSC_FALSE : PETSC_TRUE; 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2456858538eSMatthew G. Knepley } 2466858538eSMatthew G. Knepley 2476858538eSMatthew G. Knepley /*@ 24820f4b53cSBarry Smith DMGetCoordinatesLocalized - Check if the `DM` coordinates have been localized for cells 2496858538eSMatthew G. Knepley 25020f4b53cSBarry Smith Collective 2516858538eSMatthew G. Knepley 2526858538eSMatthew G. Knepley Input Parameter: 25320f4b53cSBarry Smith . dm - The `DM` 2546858538eSMatthew G. Knepley 2556858538eSMatthew G. Knepley Output Parameter: 256a4e35b19SJacob Faibussowitsch . areLocalized - `PETSC_TRUE` if localized 2576858538eSMatthew G. Knepley 2586858538eSMatthew G. Knepley Level: developer 2596858538eSMatthew G. Knepley 26020f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMSetPeriodicity()`, `DMGetCoordinatesLocalizedLocal()` 2616858538eSMatthew G. Knepley @*/ 262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalized(DM dm, PetscBool *areLocalized) 263d71ae5a4SJacob Faibussowitsch { 2646858538eSMatthew G. Knepley PetscBool localized; 2656858538eSMatthew G. Knepley 2666858538eSMatthew G. Knepley PetscFunctionBegin; 2676858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2684f572ea9SToby Isaac PetscAssertPointer(areLocalized, 2); 2696858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalizedLocal(dm, &localized)); 270462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&localized, areLocalized, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm))); 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2726858538eSMatthew G. Knepley } 2736858538eSMatthew G. Knepley 2746858538eSMatthew G. Knepley /*@ 275a4b8866cSMatthew G. Knepley DMGetSparseLocalize - Check if the `DM` coordinates should be localized only for cells near the periodic boundary. 276a4b8866cSMatthew G. Knepley 277a4b8866cSMatthew G. Knepley Not collective 278a4b8866cSMatthew G. Knepley 279a4b8866cSMatthew G. Knepley Input Parameter: 280a4b8866cSMatthew G. Knepley . dm - The `DM` 281a4b8866cSMatthew G. Knepley 282a4b8866cSMatthew G. Knepley Output Parameter: 283d7c1f440SPierre Jolivet . sparse - `PETSC_TRUE` if only cells near the periodic boundary are localized 284a4b8866cSMatthew G. Knepley 285a4b8866cSMatthew G. Knepley Level: intermediate 286a4b8866cSMatthew G. Knepley 287a4b8866cSMatthew G. Knepley .seealso: `DMSetSparseLocalize()`, `DMLocalizeCoordinates()`, `DMSetPeriodicity()` 288a4b8866cSMatthew G. Knepley @*/ 289a4b8866cSMatthew G. Knepley PetscErrorCode DMGetSparseLocalize(DM dm, PetscBool *sparse) 290a4b8866cSMatthew G. Knepley { 291a4b8866cSMatthew G. Knepley PetscFunctionBegin; 292a4b8866cSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 293a4b8866cSMatthew G. Knepley PetscAssertPointer(sparse, 2); 294a4b8866cSMatthew G. Knepley *sparse = dm->sparseLocalize; 295a4b8866cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 296a4b8866cSMatthew G. Knepley } 297a4b8866cSMatthew G. Knepley 298a4b8866cSMatthew G. Knepley /*@ 299a4b8866cSMatthew G. Knepley DMSetSparseLocalize - Set the flag indicating that `DM` coordinates should be localized only for cells near the periodic boundary. 300a4b8866cSMatthew G. Knepley 301a4b8866cSMatthew G. Knepley Logically collective 302a4b8866cSMatthew G. Knepley 303a4b8866cSMatthew G. Knepley Input Parameters: 304a4b8866cSMatthew G. Knepley + dm - The `DM` 305d7c1f440SPierre Jolivet - sparse - `PETSC_TRUE` if only cells near the periodic boundary are localized 306a4b8866cSMatthew G. Knepley 307a4b8866cSMatthew G. Knepley Level: intermediate 308a4b8866cSMatthew G. Knepley 309a4b8866cSMatthew G. Knepley .seealso: `DMGetSparseLocalize()`, `DMLocalizeCoordinates()`, `DMSetPeriodicity()` 310a4b8866cSMatthew G. Knepley @*/ 311a4b8866cSMatthew G. Knepley PetscErrorCode DMSetSparseLocalize(DM dm, PetscBool sparse) 312a4b8866cSMatthew G. Knepley { 313a4b8866cSMatthew G. Knepley PetscFunctionBegin; 314a4b8866cSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 315a4b8866cSMatthew G. Knepley PetscValidLogicalCollectiveBool(dm, sparse, 2); 316a4b8866cSMatthew G. Knepley dm->sparseLocalize = sparse; 317a4b8866cSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 318a4b8866cSMatthew G. Knepley } 319a4b8866cSMatthew G. Knepley 320a4b8866cSMatthew G. Knepley /*@ 3216858538eSMatthew G. Knepley DMLocalizeCoordinates - If a mesh is periodic, create local coordinates for cells having periodic faces 3226858538eSMatthew G. Knepley 32320f4b53cSBarry Smith Collective 3246858538eSMatthew G. Knepley 3256858538eSMatthew G. Knepley Input Parameter: 32620f4b53cSBarry Smith . dm - The `DM` 3276858538eSMatthew G. Knepley 3286858538eSMatthew G. Knepley Level: developer 3296858538eSMatthew G. Knepley 33020f4b53cSBarry Smith .seealso: `DM`, `DMSetPeriodicity()`, `DMLocalizeCoordinate()`, `DMLocalizeAddCoordinate()` 3316858538eSMatthew G. Knepley @*/ 332d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinates(DM dm) 333d71ae5a4SJacob Faibussowitsch { 3346858538eSMatthew G. Knepley DM cdm, cdgdm, cplex, plex; 3356858538eSMatthew G. Knepley PetscSection cs, csDG; 3366858538eSMatthew G. Knepley Vec coordinates, cVec; 3376858538eSMatthew G. Knepley PetscScalar *coordsDG, *anchor, *localized; 3384fb89dddSMatthew G. Knepley const PetscReal *Lstart, *L; 3391690c2aeSBarry Smith PetscInt Nc, vStart, vEnd, sStart, sEnd, newStart = PETSC_INT_MAX, newEnd = PETSC_INT_MIN, bs, coordSize; 340a4b8866cSMatthew G. Knepley PetscBool isLocalized, sparseLocalize, useDG = PETSC_FALSE, useDGGlobal; 3416858538eSMatthew G. Knepley PetscInt maxHeight = 0, h; 3426858538eSMatthew G. Knepley PetscInt *pStart = NULL, *pEnd = NULL; 3436858538eSMatthew G. Knepley MPI_Comm comm; 3446858538eSMatthew G. Knepley 3456858538eSMatthew G. Knepley PetscFunctionBegin; 3466858538eSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3474fb89dddSMatthew G. Knepley PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L)); 348a4b8866cSMatthew G. Knepley PetscCall(DMGetSparseLocalize(dm, &sparseLocalize)); 3496858538eSMatthew G. Knepley /* Cannot automatically localize without L and maxCell right now */ 3503ba16761SJacob Faibussowitsch if (!L) PetscFunctionReturn(PETSC_SUCCESS); 3516858538eSMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3526858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalized(dm, &isLocalized)); 3533ba16761SJacob Faibussowitsch if (isLocalized) PetscFunctionReturn(PETSC_SUCCESS); 3546858538eSMatthew G. Knepley 3556858538eSMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); 3566858538eSMatthew G. Knepley PetscCall(DMConvert(dm, DMPLEX, &plex)); 3576858538eSMatthew G. Knepley PetscCall(DMConvert(cdm, DMPLEX, &cplex)); 3586858538eSMatthew G. Knepley if (cplex) { 3596858538eSMatthew G. Knepley PetscCall(DMPlexGetDepthStratum(cplex, 0, &vStart, &vEnd)); 3606858538eSMatthew G. Knepley PetscCall(DMPlexGetMaxProjectionHeight(cplex, &maxHeight)); 3616858538eSMatthew G. Knepley PetscCall(DMGetWorkArray(dm, 2 * (maxHeight + 1), MPIU_INT, &pStart)); 3626858538eSMatthew G. Knepley pEnd = &pStart[maxHeight + 1]; 3636858538eSMatthew G. Knepley newStart = vStart; 3646858538eSMatthew G. Knepley newEnd = vEnd; 3656858538eSMatthew G. Knepley for (h = 0; h <= maxHeight; h++) { 3666858538eSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(cplex, h, &pStart[h], &pEnd[h])); 3676858538eSMatthew G. Knepley newStart = PetscMin(newStart, pStart[h]); 3686858538eSMatthew G. Knepley newEnd = PetscMax(newEnd, pEnd[h]); 3696858538eSMatthew G. Knepley } 3706858538eSMatthew G. Knepley } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM"); 3716858538eSMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 3726858538eSMatthew G. Knepley PetscCheck(coordinates, comm, PETSC_ERR_SUP, "Missing local coordinate vector"); 3736858538eSMatthew G. Knepley PetscCall(DMGetCoordinateSection(dm, &cs)); 3746858538eSMatthew G. Knepley PetscCall(VecGetBlockSize(coordinates, &bs)); 3756858538eSMatthew G. Knepley PetscCall(PetscSectionGetChart(cs, &sStart, &sEnd)); 3766858538eSMatthew G. Knepley 3776858538eSMatthew G. Knepley PetscCall(PetscSectionCreate(comm, &csDG)); 3786858538eSMatthew G. Knepley PetscCall(PetscSectionSetNumFields(csDG, 1)); 3796858538eSMatthew G. Knepley PetscCall(PetscSectionGetFieldComponents(cs, 0, &Nc)); 3806858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldComponents(csDG, 0, Nc)); 3816858538eSMatthew G. Knepley PetscCall(PetscSectionSetChart(csDG, newStart, newEnd)); 3826858538eSMatthew G. Knepley PetscCheck(bs == Nc, comm, PETSC_ERR_ARG_INCOMP, "Coordinate block size %" PetscInt_FMT " != %" PetscInt_FMT " number of components", bs, Nc); 3836858538eSMatthew G. Knepley 3846858538eSMatthew G. Knepley PetscCall(DMGetWorkArray(dm, 2 * Nc, MPIU_SCALAR, &anchor)); 3856858538eSMatthew G. Knepley localized = &anchor[Nc]; 3866858538eSMatthew G. Knepley for (h = 0; h <= maxHeight; h++) { 3876858538eSMatthew G. Knepley PetscInt cStart = pStart[h], cEnd = pEnd[h], c; 3886858538eSMatthew G. Knepley 3896858538eSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 3906858538eSMatthew G. Knepley PetscScalar *cellCoords = NULL; 3916858538eSMatthew G. Knepley DMPolytopeType ct; 3926858538eSMatthew G. Knepley PetscInt dof, d, p; 3936858538eSMatthew G. Knepley 3946858538eSMatthew G. Knepley PetscCall(DMPlexGetCellType(plex, c, &ct)); 3956858538eSMatthew G. Knepley if (ct == DM_POLYTOPE_FV_GHOST) continue; 3966858538eSMatthew G. Knepley PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 3976858538eSMatthew G. Knepley PetscCheck(!(dof % Nc), comm, PETSC_ERR_ARG_INCOMP, "Coordinate size on cell %" PetscInt_FMT " closure %" PetscInt_FMT " not divisible by %" PetscInt_FMT " number of components", c, dof, Nc); 3986858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[d]; 3996858538eSMatthew G. Knepley for (p = 0; p < dof / Nc; ++p) { 4006858538eSMatthew G. Knepley PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p * Nc], localized)); 4019371c9d4SSatish Balay for (d = 0; d < Nc; ++d) 4029371c9d4SSatish Balay if (cellCoords[p * Nc + d] != localized[d]) break; 4036858538eSMatthew G. Knepley if (d < Nc) break; 4046858538eSMatthew G. Knepley } 4056858538eSMatthew G. Knepley if (p < dof / Nc) useDG = PETSC_TRUE; 4066858538eSMatthew G. Knepley if (p < dof / Nc || !sparseLocalize) { 4076858538eSMatthew G. Knepley PetscCall(PetscSectionSetDof(csDG, c, dof)); 4086858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(csDG, c, 0, dof)); 4096858538eSMatthew G. Knepley } 4106858538eSMatthew G. Knepley PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 4116858538eSMatthew G. Knepley } 4126858538eSMatthew G. Knepley } 413462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&useDG, &useDGGlobal, 1, MPIU_BOOL, MPI_LOR, comm)); 4146858538eSMatthew G. Knepley if (!useDGGlobal) goto end; 4156858538eSMatthew G. Knepley 4166858538eSMatthew G. Knepley PetscCall(PetscSectionSetUp(csDG)); 4176858538eSMatthew G. Knepley PetscCall(PetscSectionGetStorageSize(csDG, &coordSize)); 4186858538eSMatthew G. Knepley PetscCall(VecCreate(PETSC_COMM_SELF, &cVec)); 4196858538eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)cVec, "coordinates")); 4206858538eSMatthew G. Knepley PetscCall(VecSetBlockSize(cVec, bs)); 4216858538eSMatthew G. Knepley PetscCall(VecSetSizes(cVec, coordSize, PETSC_DETERMINE)); 4226858538eSMatthew G. Knepley PetscCall(VecSetType(cVec, VECSTANDARD)); 4236858538eSMatthew G. Knepley PetscCall(VecGetArray(cVec, &coordsDG)); 4246858538eSMatthew G. Knepley for (h = 0; h <= maxHeight; h++) { 4256858538eSMatthew G. Knepley PetscInt cStart = pStart[h], cEnd = pEnd[h], c; 4266858538eSMatthew G. Knepley 4276858538eSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 4286858538eSMatthew G. Knepley PetscScalar *cellCoords = NULL; 4296858538eSMatthew G. Knepley PetscInt p = 0, q, dof, cdof, d, offDG; 4306858538eSMatthew G. Knepley 4316858538eSMatthew G. Knepley PetscCall(PetscSectionGetDof(csDG, c, &cdof)); 4326858538eSMatthew G. Knepley if (!cdof) continue; 4336858538eSMatthew G. Knepley PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 4346858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(csDG, c, &offDG)); 4357c48043bSMatthew G. Knepley // TODO The coordinates are set in closure order, which might not be the tensor order 4366858538eSMatthew G. Knepley for (q = 0; q < dof / Nc; ++q) { 4376858538eSMatthew G. Knepley // Select a trial anchor 4386858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[q * Nc + d]; 4396858538eSMatthew G. Knepley for (p = 0; p < dof / Nc; ++p) { 4406858538eSMatthew G. Knepley PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p * Nc], &coordsDG[offDG + p * Nc])); 4414fb89dddSMatthew G. Knepley // We need the cell to fit into the torus [lower, lower+L) 4426858538eSMatthew G. Knepley for (d = 0; d < Nc; ++d) 4434fb89dddSMatthew G. Knepley if (L[d] > 0. && ((PetscRealPart(coordsDG[offDG + p * Nc + d]) < (Lstart ? Lstart[d] : 0.)) || (PetscRealPart(coordsDG[offDG + p * Nc + d]) > (Lstart ? Lstart[d] : 0.) + L[d]))) break; 4446858538eSMatthew G. Knepley if (d < Nc) break; 4456858538eSMatthew G. Knepley } 4466858538eSMatthew G. Knepley if (p == dof / Nc) break; 4476858538eSMatthew G. Knepley } 4484fb89dddSMatthew G. Knepley PetscCheck(p == dof / Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " does not fit into the torus %s[0, L]", c, Lstart ? "Lstart + " : ""); 4496858538eSMatthew G. Knepley PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords)); 4506858538eSMatthew G. Knepley } 4516858538eSMatthew G. Knepley } 4526858538eSMatthew G. Knepley PetscCall(VecRestoreArray(cVec, &coordsDG)); 45399acd26cSksagiyam PetscUseTypeMethod(dm, createcellcoordinatedm, &cdgdm); 4546858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateDM(dm, cdgdm)); 4554c712d99Sksagiyam PetscCall(DMDestroy(&cdgdm)); 4567c48043bSMatthew G. Knepley // Convert the discretization 4577c48043bSMatthew G. Knepley { 4584c712d99Sksagiyam PetscFE fe; 4597c48043bSMatthew G. Knepley PetscClassId id; 4607c48043bSMatthew G. Knepley 4617c48043bSMatthew G. Knepley PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe)); 4627c48043bSMatthew G. Knepley PetscCall(PetscObjectGetClassId((PetscObject)fe, &id)); 4637c48043bSMatthew G. Knepley if (id == PETSCFE_CLASSID) { 4644c712d99Sksagiyam PetscSpace P; 4654c712d99Sksagiyam PetscInt degree; 4664c712d99Sksagiyam 4677c48043bSMatthew G. Knepley PetscCall(PetscFEGetBasisSpace(fe, &P)); 4684c712d99Sksagiyam PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 469*e65c294aSksagiyam PetscCall(DMPlexCreateCoordinateSpace(dm, degree, PETSC_TRUE, PETSC_FALSE)); 4707c48043bSMatthew G. Knepley } 4717c48043bSMatthew G. Knepley } 4724c712d99Sksagiyam PetscCall(DMSetCellCoordinateSection(dm, PETSC_DETERMINE, csDG)); 4734c712d99Sksagiyam PetscCall(DMSetCellCoordinatesLocal(dm, cVec)); 4744c712d99Sksagiyam PetscCall(VecDestroy(&cVec)); 4756858538eSMatthew G. Knepley 4766858538eSMatthew G. Knepley end: 4776858538eSMatthew G. Knepley PetscCall(DMRestoreWorkArray(dm, 2 * bs, MPIU_SCALAR, &anchor)); 4786858538eSMatthew G. Knepley PetscCall(DMRestoreWorkArray(dm, 2 * (maxHeight + 1), MPIU_INT, &pStart)); 4796858538eSMatthew G. Knepley PetscCall(PetscSectionDestroy(&csDG)); 4806858538eSMatthew G. Knepley PetscCall(DMDestroy(&plex)); 4816858538eSMatthew G. Knepley PetscCall(DMDestroy(&cplex)); 4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4836858538eSMatthew G. Knepley } 484