xref: /petsc/src/dm/interface/dmperiodicity.c (revision 1690c2ae071c7584458d4e437df7b47bc4686b3c)
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 
86858538eSMatthew G. Knepley   Input Parameter:
920f4b53cSBarry Smith . dm - The `DM` object
106858538eSMatthew G. Knepley 
116858538eSMatthew G. Knepley   Output Parameters:
126858538eSMatthew 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
1320f4b53cSBarry Smith . Lstart  - If we assume the mesh is a torus, this is the start of each coordinate, or `NULL` for 0.0
146858538eSMatthew G. Knepley - L       - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0
156858538eSMatthew G. Knepley 
166858538eSMatthew G. Knepley   Level: developer
176858538eSMatthew G. Knepley 
1842747ad1SJacob Faibussowitsch .seealso: `DM`
196858538eSMatthew G. Knepley @*/
205d83a8b1SBarry Smith PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal *maxCell[], const PetscReal *Lstart[], const PetscReal *L[])
21d71ae5a4SJacob Faibussowitsch {
226858538eSMatthew G. Knepley   PetscFunctionBegin;
236858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
246858538eSMatthew G. Knepley   if (maxCell) *maxCell = dm->maxCell;
254fb89dddSMatthew G. Knepley   if (Lstart) *Lstart = dm->Lstart;
264fb89dddSMatthew G. Knepley   if (L) *L = dm->L;
273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
286858538eSMatthew G. Knepley }
296858538eSMatthew G. Knepley 
30cc4c1da9SBarry Smith /*@
316858538eSMatthew G. Knepley   DMSetPeriodicity - Set the description of mesh periodicity
326858538eSMatthew G. Knepley 
336858538eSMatthew G. Knepley   Input Parameters:
3420f4b53cSBarry Smith + dm      - The `DM` object
3520f4b53cSBarry 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.
3620f4b53cSBarry Smith . Lstart  - If we assume the mesh is a torus, this is the start of each coordinate, or `NULL` for 0.0
376858538eSMatthew G. Knepley - L       - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0
386858538eSMatthew G. Knepley 
396858538eSMatthew G. Knepley   Level: developer
406858538eSMatthew G. Knepley 
4120f4b53cSBarry Smith .seealso: `DM`, `DMGetPeriodicity()`
426858538eSMatthew G. Knepley @*/
43d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal Lstart[], const PetscReal L[])
44d71ae5a4SJacob Faibussowitsch {
456858538eSMatthew G. Knepley   PetscInt dim, d;
466858538eSMatthew G. Knepley 
476858538eSMatthew G. Knepley   PetscFunctionBegin;
486858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
494f572ea9SToby Isaac   if (maxCell) PetscAssertPointer(maxCell, 2);
504f572ea9SToby Isaac   if (Lstart) PetscAssertPointer(Lstart, 3);
514f572ea9SToby Isaac   if (L) PetscAssertPointer(L, 4);
526858538eSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
536858538eSMatthew G. Knepley   if (maxCell) {
546858538eSMatthew G. Knepley     if (!dm->maxCell) PetscCall(PetscMalloc1(dim, &dm->maxCell));
556858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) dm->maxCell[d] = maxCell[d];
566858538eSMatthew G. Knepley   } else { /* remove maxCell information to disable automatic computation of localized vertices */
576858538eSMatthew G. Knepley     PetscCall(PetscFree(dm->maxCell));
586858538eSMatthew G. Knepley     dm->maxCell = NULL;
596858538eSMatthew G. Knepley   }
604fb89dddSMatthew G. Knepley   if (Lstart) {
614fb89dddSMatthew G. Knepley     if (!dm->Lstart) PetscCall(PetscMalloc1(dim, &dm->Lstart));
624fb89dddSMatthew G. Knepley     for (d = 0; d < dim; ++d) dm->Lstart[d] = Lstart[d];
634fb89dddSMatthew G. Knepley   } else { /* remove L information to disable automatic computation of localized vertices */
644fb89dddSMatthew G. Knepley     PetscCall(PetscFree(dm->Lstart));
654fb89dddSMatthew G. Knepley     dm->Lstart = NULL;
664fb89dddSMatthew G. Knepley   }
676858538eSMatthew G. Knepley   if (L) {
686858538eSMatthew G. Knepley     if (!dm->L) PetscCall(PetscMalloc1(dim, &dm->L));
696858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) dm->L[d] = L[d];
706858538eSMatthew G. Knepley   } else { /* remove L information to disable automatic computation of localized vertices */
716858538eSMatthew G. Knepley     PetscCall(PetscFree(dm->L));
726858538eSMatthew G. Knepley     dm->L = NULL;
736858538eSMatthew G. Knepley   }
746858538eSMatthew 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");
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
766858538eSMatthew G. Knepley }
776858538eSMatthew G. Knepley 
786858538eSMatthew G. Knepley /*@
796858538eSMatthew 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.
806858538eSMatthew G. Knepley 
816858538eSMatthew G. Knepley   Input Parameters:
8220f4b53cSBarry Smith + dm       - The `DM`
836858538eSMatthew G. Knepley . in       - The input coordinate point (dim numbers)
846858538eSMatthew G. Knepley - endpoint - Include the endpoint L_i
856858538eSMatthew G. Knepley 
866858538eSMatthew G. Knepley   Output Parameter:
870b3275a6SBarry Smith . out - The localized coordinate point (dim numbers)
886858538eSMatthew G. Knepley 
896858538eSMatthew G. Knepley   Level: developer
906858538eSMatthew G. Knepley 
9120f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()`
926858538eSMatthew G. Knepley @*/
93d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinate(DM dm, const PetscScalar in[], PetscBool endpoint, PetscScalar out[])
94d71ae5a4SJacob Faibussowitsch {
956858538eSMatthew G. Knepley   PetscInt dim, d;
966858538eSMatthew G. Knepley 
976858538eSMatthew G. Knepley   PetscFunctionBegin;
986858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &dim));
996858538eSMatthew G. Knepley   if (!dm->maxCell) {
1006858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] = in[d];
1016858538eSMatthew G. Knepley   } else {
1026858538eSMatthew G. Knepley     if (endpoint) {
1036858538eSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
1046858538eSMatthew 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)) {
1056858538eSMatthew G. Knepley           out[d] = in[d] - dm->L[d] * (PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]) - 1);
1066858538eSMatthew G. Knepley         } else {
1076858538eSMatthew G. Knepley           out[d] = in[d] - dm->L[d] * PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]);
1086858538eSMatthew G. Knepley         }
1096858538eSMatthew G. Knepley       }
1106858538eSMatthew G. Knepley     } else {
111ad540459SPierre Jolivet       for (d = 0; d < dim; ++d) out[d] = in[d] - dm->L[d] * PetscFloorReal(PetscRealPart(in[d]) / dm->L[d]);
1126858538eSMatthew G. Knepley     }
1136858538eSMatthew G. Knepley   }
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1156858538eSMatthew G. Knepley }
1166858538eSMatthew G. Knepley 
1176858538eSMatthew G. Knepley /*
1186858538eSMatthew 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.
1196858538eSMatthew G. Knepley 
1206858538eSMatthew G. Knepley   Input Parameters:
12120f4b53cSBarry Smith + dm     - The `DM`
1226858538eSMatthew G. Knepley . dim    - The spatial dimension
1236858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it
1246858538eSMatthew G. Knepley - in     - The input coordinate point (dim numbers)
1256858538eSMatthew G. Knepley 
1266858538eSMatthew G. Knepley   Output Parameter:
1270b3275a6SBarry Smith . out - The localized coordinate point (dim numbers)
1286858538eSMatthew G. Knepley 
1296858538eSMatthew G. Knepley   Level: developer
1306858538eSMatthew G. Knepley 
13120f4b53cSBarry Smith   Note:
13220f4b53cSBarry 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
1336858538eSMatthew G. Knepley 
13420f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()`
1356858538eSMatthew G. Knepley */
136d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[])
137d71ae5a4SJacob Faibussowitsch {
1386858538eSMatthew G. Knepley   PetscInt d;
1396858538eSMatthew G. Knepley 
1406858538eSMatthew G. Knepley   PetscFunctionBegin;
1416858538eSMatthew G. Knepley   if (!dm->maxCell) {
1426858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] = in[d];
1436858538eSMatthew G. Knepley   } else {
1446858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1456858538eSMatthew G. Knepley       if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) {
1466858538eSMatthew G. Knepley         out[d] = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d];
1476858538eSMatthew G. Knepley       } else {
1486858538eSMatthew G. Knepley         out[d] = in[d];
1496858538eSMatthew G. Knepley       }
1506858538eSMatthew G. Knepley     }
1516858538eSMatthew G. Knepley   }
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1536858538eSMatthew G. Knepley }
1546858538eSMatthew G. Knepley 
155d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim, const PetscReal anchor[], const PetscReal in[], PetscReal out[])
156d71ae5a4SJacob Faibussowitsch {
1576858538eSMatthew G. Knepley   PetscInt d;
1586858538eSMatthew G. Knepley 
1596858538eSMatthew G. Knepley   PetscFunctionBegin;
1606858538eSMatthew G. Knepley   if (!dm->maxCell) {
1616858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] = in[d];
1626858538eSMatthew G. Knepley   } else {
1636858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1646858538eSMatthew G. Knepley       if ((dm->L[d] > 0.0) && (PetscAbsReal(anchor[d] - in[d]) > dm->maxCell[d])) {
1656858538eSMatthew G. Knepley         out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d];
1666858538eSMatthew G. Knepley       } else {
1676858538eSMatthew G. Knepley         out[d] = in[d];
1686858538eSMatthew G. Knepley       }
1696858538eSMatthew G. Knepley     }
1706858538eSMatthew G. Knepley   }
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1726858538eSMatthew G. Knepley }
1736858538eSMatthew G. Knepley 
1746858538eSMatthew G. Knepley /*
1756858538eSMatthew 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.
1766858538eSMatthew G. Knepley 
1776858538eSMatthew G. Knepley   Input Parameters:
17820f4b53cSBarry Smith + dm     - The `DM`
1796858538eSMatthew G. Knepley . dim    - The spatial dimension
1806858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it
1816858538eSMatthew G. Knepley . in     - The input coordinate delta (dim numbers)
1826858538eSMatthew G. Knepley - out    - The input coordinate point (dim numbers)
1836858538eSMatthew G. Knepley 
1846858538eSMatthew G. Knepley   Output Parameter:
1856858538eSMatthew G. Knepley . out    - The localized coordinate in + out
1866858538eSMatthew G. Knepley 
1876858538eSMatthew G. Knepley   Level: developer
1886858538eSMatthew G. Knepley 
18920f4b53cSBarry Smith   Note:
1900b3275a6SBarry 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
1916858538eSMatthew G. Knepley 
19220f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMLocalizeCoordinate()`
1936858538eSMatthew G. Knepley */
194d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeAddCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[])
195d71ae5a4SJacob Faibussowitsch {
1966858538eSMatthew G. Knepley   PetscInt d;
1976858538eSMatthew G. Knepley 
1986858538eSMatthew G. Knepley   PetscFunctionBegin;
1996858538eSMatthew G. Knepley   if (!dm->maxCell) {
2006858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] += in[d];
2016858538eSMatthew G. Knepley   } else {
2026858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
2036858538eSMatthew G. Knepley       const PetscReal maxC = dm->maxCell[d];
2046858538eSMatthew G. Knepley 
2056858538eSMatthew G. Knepley       if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > maxC)) {
2066858538eSMatthew G. Knepley         const PetscScalar newCoord = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d];
2076858538eSMatthew G. Knepley 
2086858538eSMatthew G. Knepley         if (PetscAbsScalar(newCoord - anchor[d]) > maxC)
2096858538eSMatthew 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]));
2106858538eSMatthew G. Knepley         out[d] += newCoord;
2116858538eSMatthew G. Knepley       } else {
2126858538eSMatthew G. Knepley         out[d] += in[d];
2136858538eSMatthew G. Knepley       }
2146858538eSMatthew G. Knepley     }
2156858538eSMatthew G. Knepley   }
2163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2176858538eSMatthew G. Knepley }
2186858538eSMatthew G. Knepley 
2196858538eSMatthew G. Knepley /*@
22020f4b53cSBarry Smith   DMGetCoordinatesLocalizedLocal - Check if the `DM` coordinates have been localized for cells on this process
2216858538eSMatthew G. Knepley 
22220f4b53cSBarry Smith   Not Collective
2236858538eSMatthew G. Knepley 
2246858538eSMatthew G. Knepley   Input Parameter:
22520f4b53cSBarry Smith . dm - The `DM`
2266858538eSMatthew G. Knepley 
2276858538eSMatthew G. Knepley   Output Parameter:
228a4e35b19SJacob Faibussowitsch . areLocalized - `PETSC_TRUE` if localized
2296858538eSMatthew G. Knepley 
2306858538eSMatthew G. Knepley   Level: developer
2316858538eSMatthew G. Knepley 
23220f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMGetCoordinatesLocalized()`, `DMSetPeriodicity()`
2336858538eSMatthew G. Knepley @*/
234d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalizedLocal(DM dm, PetscBool *areLocalized)
235d71ae5a4SJacob Faibussowitsch {
2366858538eSMatthew G. Knepley   PetscFunctionBegin;
2376858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2384f572ea9SToby Isaac   PetscAssertPointer(areLocalized, 2);
2396858538eSMatthew G. Knepley   *areLocalized = dm->coordinates[1].dim < 0 ? PETSC_FALSE : PETSC_TRUE;
2403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2416858538eSMatthew G. Knepley }
2426858538eSMatthew G. Knepley 
2436858538eSMatthew G. Knepley /*@
24420f4b53cSBarry Smith   DMGetCoordinatesLocalized - Check if the `DM` coordinates have been localized for cells
2456858538eSMatthew G. Knepley 
24620f4b53cSBarry Smith   Collective
2476858538eSMatthew G. Knepley 
2486858538eSMatthew G. Knepley   Input Parameter:
24920f4b53cSBarry Smith . dm - The `DM`
2506858538eSMatthew G. Knepley 
2516858538eSMatthew G. Knepley   Output Parameter:
252a4e35b19SJacob Faibussowitsch . areLocalized - `PETSC_TRUE` if localized
2536858538eSMatthew G. Knepley 
2546858538eSMatthew G. Knepley   Level: developer
2556858538eSMatthew G. Knepley 
25620f4b53cSBarry Smith .seealso: `DM`, `DMLocalizeCoordinates()`, `DMSetPeriodicity()`, `DMGetCoordinatesLocalizedLocal()`
2576858538eSMatthew G. Knepley @*/
258d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalized(DM dm, PetscBool *areLocalized)
259d71ae5a4SJacob Faibussowitsch {
2606858538eSMatthew G. Knepley   PetscBool localized;
2616858538eSMatthew G. Knepley 
2626858538eSMatthew G. Knepley   PetscFunctionBegin;
2636858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2644f572ea9SToby Isaac   PetscAssertPointer(areLocalized, 2);
2656858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalizedLocal(dm, &localized));
2666858538eSMatthew G. Knepley   PetscCall(MPIU_Allreduce(&localized, areLocalized, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
2673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2686858538eSMatthew G. Knepley }
2696858538eSMatthew G. Knepley 
2706858538eSMatthew G. Knepley /*@
271a4b8866cSMatthew G. Knepley   DMGetSparseLocalize - Check if the `DM` coordinates should be localized only for cells near the periodic boundary.
272a4b8866cSMatthew G. Knepley 
273a4b8866cSMatthew G. Knepley   Not collective
274a4b8866cSMatthew G. Knepley 
275a4b8866cSMatthew G. Knepley   Input Parameter:
276a4b8866cSMatthew G. Knepley . dm - The `DM`
277a4b8866cSMatthew G. Knepley 
278a4b8866cSMatthew G. Knepley   Output Parameter:
279a4b8866cSMatthew G. Knepley . sparse - `PETSC_TRUE` if ony cells near the periodic boundary are localized
280a4b8866cSMatthew G. Knepley 
281a4b8866cSMatthew G. Knepley   Level: intermediate
282a4b8866cSMatthew G. Knepley 
283a4b8866cSMatthew G. Knepley .seealso: `DMSetSparseLocalize()`, `DMLocalizeCoordinates()`, `DMSetPeriodicity()`
284a4b8866cSMatthew G. Knepley @*/
285a4b8866cSMatthew G. Knepley PetscErrorCode DMGetSparseLocalize(DM dm, PetscBool *sparse)
286a4b8866cSMatthew G. Knepley {
287a4b8866cSMatthew G. Knepley   PetscFunctionBegin;
288a4b8866cSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
289a4b8866cSMatthew G. Knepley   PetscAssertPointer(sparse, 2);
290a4b8866cSMatthew G. Knepley   *sparse = dm->sparseLocalize;
291a4b8866cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
292a4b8866cSMatthew G. Knepley }
293a4b8866cSMatthew G. Knepley 
294a4b8866cSMatthew G. Knepley /*@
295a4b8866cSMatthew G. Knepley   DMSetSparseLocalize - Set the flag indicating that `DM` coordinates should be localized only for cells near the periodic boundary.
296a4b8866cSMatthew G. Knepley 
297a4b8866cSMatthew G. Knepley   Logically collective
298a4b8866cSMatthew G. Knepley 
299a4b8866cSMatthew G. Knepley   Input Parameters:
300a4b8866cSMatthew G. Knepley + dm     - The `DM`
301a4b8866cSMatthew G. Knepley - sparse - `PETSC_TRUE` if ony cells near the periodic boundary are localized
302a4b8866cSMatthew G. Knepley 
303a4b8866cSMatthew G. Knepley   Level: intermediate
304a4b8866cSMatthew G. Knepley 
305a4b8866cSMatthew G. Knepley .seealso: `DMGetSparseLocalize()`, `DMLocalizeCoordinates()`, `DMSetPeriodicity()`
306a4b8866cSMatthew G. Knepley @*/
307a4b8866cSMatthew G. Knepley PetscErrorCode DMSetSparseLocalize(DM dm, PetscBool sparse)
308a4b8866cSMatthew G. Knepley {
309a4b8866cSMatthew G. Knepley   PetscFunctionBegin;
310a4b8866cSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
311a4b8866cSMatthew G. Knepley   PetscValidLogicalCollectiveBool(dm, sparse, 2);
312a4b8866cSMatthew G. Knepley   dm->sparseLocalize = sparse;
313a4b8866cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
314a4b8866cSMatthew G. Knepley }
315a4b8866cSMatthew G. Knepley 
316a4b8866cSMatthew G. Knepley /*@
3176858538eSMatthew G. Knepley   DMLocalizeCoordinates - If a mesh is periodic, create local coordinates for cells having periodic faces
3186858538eSMatthew G. Knepley 
31920f4b53cSBarry Smith   Collective
3206858538eSMatthew G. Knepley 
3216858538eSMatthew G. Knepley   Input Parameter:
32220f4b53cSBarry Smith . dm - The `DM`
3236858538eSMatthew G. Knepley 
3246858538eSMatthew G. Knepley   Level: developer
3256858538eSMatthew G. Knepley 
32620f4b53cSBarry Smith .seealso: `DM`, `DMSetPeriodicity()`, `DMLocalizeCoordinate()`, `DMLocalizeAddCoordinate()`
3276858538eSMatthew G. Knepley @*/
328d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalizeCoordinates(DM dm)
329d71ae5a4SJacob Faibussowitsch {
3306858538eSMatthew G. Knepley   DM               cdm, cdgdm, cplex, plex;
3316858538eSMatthew G. Knepley   PetscSection     cs, csDG;
3326858538eSMatthew G. Knepley   Vec              coordinates, cVec;
3336858538eSMatthew G. Knepley   PetscScalar     *coordsDG, *anchor, *localized;
3344fb89dddSMatthew G. Knepley   const PetscReal *Lstart, *L;
335*1690c2aeSBarry Smith   PetscInt         Nc, vStart, vEnd, sStart, sEnd, newStart = PETSC_INT_MAX, newEnd = PETSC_INT_MIN, bs, coordSize;
336a4b8866cSMatthew G. Knepley   PetscBool        isLocalized, sparseLocalize, useDG = PETSC_FALSE, useDGGlobal;
3376858538eSMatthew G. Knepley   PetscInt         maxHeight = 0, h;
3386858538eSMatthew G. Knepley   PetscInt        *pStart = NULL, *pEnd = NULL;
3396858538eSMatthew G. Knepley   MPI_Comm         comm;
3406858538eSMatthew G. Knepley 
3416858538eSMatthew G. Knepley   PetscFunctionBegin;
3426858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3434fb89dddSMatthew G. Knepley   PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
344a4b8866cSMatthew G. Knepley   PetscCall(DMGetSparseLocalize(dm, &sparseLocalize));
3456858538eSMatthew G. Knepley   /* Cannot automatically localize without L and maxCell right now */
3463ba16761SJacob Faibussowitsch   if (!L) PetscFunctionReturn(PETSC_SUCCESS);
3476858538eSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3486858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalized(dm, &isLocalized));
3493ba16761SJacob Faibussowitsch   if (isLocalized) PetscFunctionReturn(PETSC_SUCCESS);
3506858538eSMatthew G. Knepley 
3516858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
3526858538eSMatthew G. Knepley   PetscCall(DMConvert(dm, DMPLEX, &plex));
3536858538eSMatthew G. Knepley   PetscCall(DMConvert(cdm, DMPLEX, &cplex));
3546858538eSMatthew G. Knepley   if (cplex) {
3556858538eSMatthew G. Knepley     PetscCall(DMPlexGetDepthStratum(cplex, 0, &vStart, &vEnd));
3566858538eSMatthew G. Knepley     PetscCall(DMPlexGetMaxProjectionHeight(cplex, &maxHeight));
3576858538eSMatthew G. Knepley     PetscCall(DMGetWorkArray(dm, 2 * (maxHeight + 1), MPIU_INT, &pStart));
3586858538eSMatthew G. Knepley     pEnd     = &pStart[maxHeight + 1];
3596858538eSMatthew G. Knepley     newStart = vStart;
3606858538eSMatthew G. Knepley     newEnd   = vEnd;
3616858538eSMatthew G. Knepley     for (h = 0; h <= maxHeight; h++) {
3626858538eSMatthew G. Knepley       PetscCall(DMPlexGetHeightStratum(cplex, h, &pStart[h], &pEnd[h]));
3636858538eSMatthew G. Knepley       newStart = PetscMin(newStart, pStart[h]);
3646858538eSMatthew G. Knepley       newEnd   = PetscMax(newEnd, pEnd[h]);
3656858538eSMatthew G. Knepley     }
3666858538eSMatthew G. Knepley   } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM");
3676858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
3686858538eSMatthew G. Knepley   PetscCheck(coordinates, comm, PETSC_ERR_SUP, "Missing local coordinate vector");
3696858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateSection(dm, &cs));
3706858538eSMatthew G. Knepley   PetscCall(VecGetBlockSize(coordinates, &bs));
3716858538eSMatthew G. Knepley   PetscCall(PetscSectionGetChart(cs, &sStart, &sEnd));
3726858538eSMatthew G. Knepley 
3736858538eSMatthew G. Knepley   PetscCall(PetscSectionCreate(comm, &csDG));
3746858538eSMatthew G. Knepley   PetscCall(PetscSectionSetNumFields(csDG, 1));
3756858538eSMatthew G. Knepley   PetscCall(PetscSectionGetFieldComponents(cs, 0, &Nc));
3766858538eSMatthew G. Knepley   PetscCall(PetscSectionSetFieldComponents(csDG, 0, Nc));
3776858538eSMatthew G. Knepley   PetscCall(PetscSectionSetChart(csDG, newStart, newEnd));
3786858538eSMatthew G. Knepley   PetscCheck(bs == Nc, comm, PETSC_ERR_ARG_INCOMP, "Coordinate block size %" PetscInt_FMT " != %" PetscInt_FMT " number of components", bs, Nc);
3796858538eSMatthew G. Knepley 
3806858538eSMatthew G. Knepley   PetscCall(DMGetWorkArray(dm, 2 * Nc, MPIU_SCALAR, &anchor));
3816858538eSMatthew G. Knepley   localized = &anchor[Nc];
3826858538eSMatthew G. Knepley   for (h = 0; h <= maxHeight; h++) {
3836858538eSMatthew G. Knepley     PetscInt cStart = pStart[h], cEnd = pEnd[h], c;
3846858538eSMatthew G. Knepley 
3856858538eSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
3866858538eSMatthew G. Knepley       PetscScalar   *cellCoords = NULL;
3876858538eSMatthew G. Knepley       DMPolytopeType ct;
3886858538eSMatthew G. Knepley       PetscInt       dof, d, p;
3896858538eSMatthew G. Knepley 
3906858538eSMatthew G. Knepley       PetscCall(DMPlexGetCellType(plex, c, &ct));
3916858538eSMatthew G. Knepley       if (ct == DM_POLYTOPE_FV_GHOST) continue;
3926858538eSMatthew G. Knepley       PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
3936858538eSMatthew 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);
3946858538eSMatthew G. Knepley       for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[d];
3956858538eSMatthew G. Knepley       for (p = 0; p < dof / Nc; ++p) {
3966858538eSMatthew G. Knepley         PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p * Nc], localized));
3979371c9d4SSatish Balay         for (d = 0; d < Nc; ++d)
3989371c9d4SSatish Balay           if (cellCoords[p * Nc + d] != localized[d]) break;
3996858538eSMatthew G. Knepley         if (d < Nc) break;
4006858538eSMatthew G. Knepley       }
4016858538eSMatthew G. Knepley       if (p < dof / Nc) useDG = PETSC_TRUE;
4026858538eSMatthew G. Knepley       if (p < dof / Nc || !sparseLocalize) {
4036858538eSMatthew G. Knepley         PetscCall(PetscSectionSetDof(csDG, c, dof));
4046858538eSMatthew G. Knepley         PetscCall(PetscSectionSetFieldDof(csDG, c, 0, dof));
4056858538eSMatthew G. Knepley       }
4066858538eSMatthew G. Knepley       PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
4076858538eSMatthew G. Knepley     }
4086858538eSMatthew G. Knepley   }
409712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&useDG, &useDGGlobal, 1, MPIU_BOOL, MPI_LOR, comm));
4106858538eSMatthew G. Knepley   if (!useDGGlobal) goto end;
4116858538eSMatthew G. Knepley 
4126858538eSMatthew G. Knepley   PetscCall(PetscSectionSetUp(csDG));
4136858538eSMatthew G. Knepley   PetscCall(PetscSectionGetStorageSize(csDG, &coordSize));
4146858538eSMatthew G. Knepley   PetscCall(VecCreate(PETSC_COMM_SELF, &cVec));
4156858538eSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)cVec, "coordinates"));
4166858538eSMatthew G. Knepley   PetscCall(VecSetBlockSize(cVec, bs));
4176858538eSMatthew G. Knepley   PetscCall(VecSetSizes(cVec, coordSize, PETSC_DETERMINE));
4186858538eSMatthew G. Knepley   PetscCall(VecSetType(cVec, VECSTANDARD));
4196858538eSMatthew G. Knepley   PetscCall(VecGetArray(cVec, &coordsDG));
4206858538eSMatthew G. Knepley   for (h = 0; h <= maxHeight; h++) {
4216858538eSMatthew G. Knepley     PetscInt cStart = pStart[h], cEnd = pEnd[h], c;
4226858538eSMatthew G. Knepley 
4236858538eSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
4246858538eSMatthew G. Knepley       PetscScalar *cellCoords = NULL;
4256858538eSMatthew G. Knepley       PetscInt     p          = 0, q, dof, cdof, d, offDG;
4266858538eSMatthew G. Knepley 
4276858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(csDG, c, &cdof));
4286858538eSMatthew G. Knepley       if (!cdof) continue;
4296858538eSMatthew G. Knepley       PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
4306858538eSMatthew G. Knepley       PetscCall(PetscSectionGetOffset(csDG, c, &offDG));
4317c48043bSMatthew G. Knepley       // TODO The coordinates are set in closure order, which might not be the tensor order
4326858538eSMatthew G. Knepley       for (q = 0; q < dof / Nc; ++q) {
4336858538eSMatthew G. Knepley         // Select a trial anchor
4346858538eSMatthew G. Knepley         for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[q * Nc + d];
4356858538eSMatthew G. Knepley         for (p = 0; p < dof / Nc; ++p) {
4366858538eSMatthew G. Knepley           PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p * Nc], &coordsDG[offDG + p * Nc]));
4374fb89dddSMatthew G. Knepley           // We need the cell to fit into the torus [lower, lower+L)
4386858538eSMatthew G. Knepley           for (d = 0; d < Nc; ++d)
4394fb89dddSMatthew 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;
4406858538eSMatthew G. Knepley           if (d < Nc) break;
4416858538eSMatthew G. Knepley         }
4426858538eSMatthew G. Knepley         if (p == dof / Nc) break;
4436858538eSMatthew G. Knepley       }
4444fb89dddSMatthew 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 + " : "");
4456858538eSMatthew G. Knepley       PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
4466858538eSMatthew G. Knepley     }
4476858538eSMatthew G. Knepley   }
4486858538eSMatthew G. Knepley   PetscCall(VecRestoreArray(cVec, &coordsDG));
4496858538eSMatthew G. Knepley   PetscCall(DMClone(cdm, &cdgdm));
4506858538eSMatthew G. Knepley   PetscCall(DMSetCellCoordinateDM(dm, cdgdm));
4516858538eSMatthew G. Knepley   PetscCall(DMSetCellCoordinateSection(dm, PETSC_DETERMINE, csDG));
4526858538eSMatthew G. Knepley   PetscCall(DMSetCellCoordinatesLocal(dm, cVec));
4536858538eSMatthew G. Knepley   PetscCall(VecDestroy(&cVec));
4547c48043bSMatthew G. Knepley   // Convert the discretization
4557c48043bSMatthew G. Knepley   {
4567c48043bSMatthew G. Knepley     PetscFE         fe, dgfe;
4577c48043bSMatthew G. Knepley     PetscSpace      P;
4587c48043bSMatthew G. Knepley     PetscDualSpace  Q, dgQ;
4597c48043bSMatthew G. Knepley     PetscQuadrature q, fq;
4607c48043bSMatthew G. Knepley     PetscClassId    id;
4617c48043bSMatthew G. Knepley 
4627c48043bSMatthew G. Knepley     PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
4637c48043bSMatthew G. Knepley     PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
4647c48043bSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
4657c48043bSMatthew G. Knepley       PetscCall(PetscFEGetBasisSpace(fe, &P));
4667c48043bSMatthew G. Knepley       PetscCall(PetscObjectReference((PetscObject)P));
4677c48043bSMatthew G. Knepley       PetscCall(PetscFEGetDualSpace(fe, &Q));
4687c48043bSMatthew G. Knepley       PetscCall(PetscDualSpaceDuplicate(Q, &dgQ));
4697c48043bSMatthew G. Knepley       PetscCall(PetscDualSpaceLagrangeSetContinuity(dgQ, PETSC_FALSE));
4707c48043bSMatthew G. Knepley       PetscCall(PetscDualSpaceSetUp(dgQ));
4717c48043bSMatthew G. Knepley       PetscCall(PetscFEGetQuadrature(fe, &q));
4727c48043bSMatthew G. Knepley       PetscCall(PetscObjectReference((PetscObject)q));
4737c48043bSMatthew G. Knepley       PetscCall(PetscFEGetFaceQuadrature(fe, &fq));
4747c48043bSMatthew G. Knepley       PetscCall(PetscObjectReference((PetscObject)fq));
4757c48043bSMatthew G. Knepley       PetscCall(PetscFECreateFromSpaces(P, dgQ, q, fq, &dgfe));
4767c48043bSMatthew G. Knepley       PetscCall(DMSetField(cdgdm, 0, NULL, (PetscObject)dgfe));
4777c48043bSMatthew G. Knepley       PetscCall(PetscFEDestroy(&dgfe));
4787c48043bSMatthew G. Knepley       PetscCall(DMCreateDS(cdgdm));
4797c48043bSMatthew G. Knepley     }
4807c48043bSMatthew G. Knepley   }
4816858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cdgdm));
4826858538eSMatthew G. Knepley 
4836858538eSMatthew G. Knepley end:
4846858538eSMatthew G. Knepley   PetscCall(DMRestoreWorkArray(dm, 2 * bs, MPIU_SCALAR, &anchor));
4856858538eSMatthew G. Knepley   PetscCall(DMRestoreWorkArray(dm, 2 * (maxHeight + 1), MPIU_INT, &pStart));
4866858538eSMatthew G. Knepley   PetscCall(PetscSectionDestroy(&csDG));
4876858538eSMatthew G. Knepley   PetscCall(DMDestroy(&plex));
4886858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cplex));
4893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4906858538eSMatthew G. Knepley }
491