xref: /petsc/src/dm/interface/dmperiodicity.c (revision 7c48043b109ee9f116d1c17e77d6369f7a4a2b2e)
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:
96858538eSMatthew G. Knepley . 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
136858538eSMatthew G. Knepley - L       - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0
146858538eSMatthew G. Knepley 
156858538eSMatthew G. Knepley   Level: developer
166858538eSMatthew G. Knepley 
176858538eSMatthew G. Knepley .seealso: `DMGetPeriodicity()`
186858538eSMatthew G. Knepley @*/
196858538eSMatthew G. Knepley PetscErrorCode DMGetPeriodicity(DM dm, const PetscReal **maxCell, const PetscReal **L)
206858538eSMatthew G. Knepley {
216858538eSMatthew G. Knepley   PetscFunctionBegin;
226858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
236858538eSMatthew G. Knepley   if (L)       *L       = dm->L;
246858538eSMatthew G. Knepley   if (maxCell) *maxCell = dm->maxCell;
256858538eSMatthew G. Knepley   PetscFunctionReturn(0);
266858538eSMatthew G. Knepley }
276858538eSMatthew G. Knepley 
286858538eSMatthew G. Knepley /*@C
296858538eSMatthew G. Knepley   DMSetPeriodicity - Set the description of mesh periodicity
306858538eSMatthew G. Knepley 
316858538eSMatthew G. Knepley   Input Parameters:
326858538eSMatthew G. Knepley + dm      - The DM object
336858538eSMatthew 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. Pass NULL to remove such information.
346858538eSMatthew G. Knepley - L       - If we assume the mesh is a torus, this is the length of each coordinate, otherwise it is < 0.0
356858538eSMatthew G. Knepley 
366858538eSMatthew G. Knepley   Level: developer
376858538eSMatthew G. Knepley 
386858538eSMatthew G. Knepley .seealso: `DMGetPeriodicity()`
396858538eSMatthew G. Knepley @*/
406858538eSMatthew G. Knepley PetscErrorCode DMSetPeriodicity(DM dm, const PetscReal maxCell[], const PetscReal L[])
416858538eSMatthew G. Knepley {
426858538eSMatthew G. Knepley   PetscInt dim, d;
436858538eSMatthew G. Knepley 
446858538eSMatthew G. Knepley   PetscFunctionBegin;
456858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
466858538eSMatthew G. Knepley   if (maxCell) {PetscValidRealPointer(maxCell,2);}
476858538eSMatthew G. Knepley   if (L)       {PetscValidRealPointer(L,3);}
486858538eSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
496858538eSMatthew G. Knepley   if (maxCell) {
506858538eSMatthew G. Knepley     if (!dm->maxCell) PetscCall(PetscMalloc1(dim, &dm->maxCell));
516858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) dm->maxCell[d] = maxCell[d];
526858538eSMatthew G. Knepley   } else { /* remove maxCell information to disable automatic computation of localized vertices */
536858538eSMatthew G. Knepley     PetscCall(PetscFree(dm->maxCell));
546858538eSMatthew G. Knepley     dm->maxCell = NULL;
556858538eSMatthew G. Knepley   }
566858538eSMatthew G. Knepley   if (L) {
576858538eSMatthew G. Knepley     if (!dm->L) PetscCall(PetscMalloc1(dim, &dm->L));
586858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) dm->L[d] = L[d];
596858538eSMatthew G. Knepley   } else { /* remove L information to disable automatic computation of localized vertices */
606858538eSMatthew G. Knepley     PetscCall(PetscFree(dm->L));
616858538eSMatthew G. Knepley     dm->L = NULL;
626858538eSMatthew G. Knepley   }
636858538eSMatthew 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");
646858538eSMatthew G. Knepley   PetscFunctionReturn(0);
656858538eSMatthew G. Knepley }
666858538eSMatthew G. Knepley 
676858538eSMatthew G. Knepley /*@
686858538eSMatthew 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.
696858538eSMatthew G. Knepley 
706858538eSMatthew G. Knepley   Input Parameters:
716858538eSMatthew G. Knepley + dm     - The DM
726858538eSMatthew G. Knepley . in     - The input coordinate point (dim numbers)
736858538eSMatthew G. Knepley - endpoint - Include the endpoint L_i
746858538eSMatthew G. Knepley 
756858538eSMatthew G. Knepley   Output Parameter:
766858538eSMatthew G. Knepley . out - The localized coordinate point
776858538eSMatthew G. Knepley 
786858538eSMatthew G. Knepley   Level: developer
796858538eSMatthew G. Knepley 
806858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()`
816858538eSMatthew G. Knepley @*/
826858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinate(DM dm, const PetscScalar in[], PetscBool endpoint, PetscScalar out[])
836858538eSMatthew G. Knepley {
846858538eSMatthew G. Knepley   PetscInt       dim, d;
856858538eSMatthew G. Knepley 
866858538eSMatthew G. Knepley   PetscFunctionBegin;
876858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &dim));
886858538eSMatthew G. Knepley   if (!dm->maxCell) {
896858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] = in[d];
906858538eSMatthew G. Knepley   } else {
916858538eSMatthew G. Knepley     if (endpoint) {
926858538eSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
936858538eSMatthew 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)) {
946858538eSMatthew G. Knepley           out[d] = in[d] - dm->L[d]*(PetscFloorReal(PetscRealPart(in[d])/dm->L[d]) - 1);
956858538eSMatthew G. Knepley         } else {
966858538eSMatthew G. Knepley           out[d] = in[d] - dm->L[d]*PetscFloorReal(PetscRealPart(in[d])/dm->L[d]);
976858538eSMatthew G. Knepley         }
986858538eSMatthew G. Knepley       }
996858538eSMatthew G. Knepley     } else {
1006858538eSMatthew G. Knepley       for (d = 0; d < dim; ++d) {
1016858538eSMatthew G. Knepley         out[d] = in[d] - dm->L[d]*PetscFloorReal(PetscRealPart(in[d])/dm->L[d]);
1026858538eSMatthew G. Knepley       }
1036858538eSMatthew G. Knepley     }
1046858538eSMatthew G. Knepley   }
1056858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1066858538eSMatthew G. Knepley }
1076858538eSMatthew G. Knepley 
1086858538eSMatthew G. Knepley /*
1096858538eSMatthew 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.
1106858538eSMatthew G. Knepley 
1116858538eSMatthew G. Knepley   Input Parameters:
1126858538eSMatthew G. Knepley + dm     - The DM
1136858538eSMatthew G. Knepley . dim    - The spatial dimension
1146858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it
1156858538eSMatthew G. Knepley - in     - The input coordinate point (dim numbers)
1166858538eSMatthew G. Knepley 
1176858538eSMatthew G. Knepley   Output Parameter:
1186858538eSMatthew G. Knepley . out - The localized coordinate point
1196858538eSMatthew G. Knepley 
1206858538eSMatthew G. Knepley   Level: developer
1216858538eSMatthew G. Knepley 
1226858538eSMatthew G. Knepley   Note: 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
1236858538eSMatthew G. Knepley 
1246858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMLocalizeAddCoordinate()`
1256858538eSMatthew G. Knepley */
1266858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[])
1276858538eSMatthew G. Knepley {
1286858538eSMatthew G. Knepley   PetscInt d;
1296858538eSMatthew G. Knepley 
1306858538eSMatthew G. Knepley   PetscFunctionBegin;
1316858538eSMatthew G. Knepley   if (!dm->maxCell) {
1326858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] = in[d];
1336858538eSMatthew G. Knepley   } else {
1346858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1356858538eSMatthew G. Knepley       if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > dm->maxCell[d])) {
1366858538eSMatthew G. Knepley         out[d] = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d];
1376858538eSMatthew G. Knepley       } else {
1386858538eSMatthew G. Knepley         out[d] = in[d];
1396858538eSMatthew G. Knepley       }
1406858538eSMatthew G. Knepley     }
1416858538eSMatthew G. Knepley   }
1426858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1436858538eSMatthew G. Knepley }
1446858538eSMatthew G. Knepley 
1456858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinateReal_Internal(DM dm, PetscInt dim, const PetscReal anchor[], const PetscReal in[], PetscReal out[])
1466858538eSMatthew G. Knepley {
1476858538eSMatthew G. Knepley   PetscInt d;
1486858538eSMatthew G. Knepley 
1496858538eSMatthew G. Knepley   PetscFunctionBegin;
1506858538eSMatthew G. Knepley   if (!dm->maxCell) {
1516858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] = in[d];
1526858538eSMatthew G. Knepley   } else {
1536858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1546858538eSMatthew G. Knepley       if ((dm->L[d] > 0.0) && (PetscAbsReal(anchor[d] - in[d]) > dm->maxCell[d])) {
1556858538eSMatthew G. Knepley         out[d] = anchor[d] > in[d] ? dm->L[d] + in[d] : in[d] - dm->L[d];
1566858538eSMatthew G. Knepley       } else {
1576858538eSMatthew G. Knepley         out[d] = in[d];
1586858538eSMatthew G. Knepley       }
1596858538eSMatthew G. Knepley     }
1606858538eSMatthew G. Knepley   }
1616858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1626858538eSMatthew G. Knepley }
1636858538eSMatthew G. Knepley 
1646858538eSMatthew G. Knepley /*
1656858538eSMatthew 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.
1666858538eSMatthew G. Knepley 
1676858538eSMatthew G. Knepley   Input Parameters:
1686858538eSMatthew G. Knepley + dm     - The DM
1696858538eSMatthew G. Knepley . dim    - The spatial dimension
1706858538eSMatthew G. Knepley . anchor - The anchor point, the input point can be no more than maxCell away from it
1716858538eSMatthew G. Knepley . in     - The input coordinate delta (dim numbers)
1726858538eSMatthew G. Knepley - out    - The input coordinate point (dim numbers)
1736858538eSMatthew G. Knepley 
1746858538eSMatthew G. Knepley   Output Parameter:
1756858538eSMatthew G. Knepley . out    - The localized coordinate in + out
1766858538eSMatthew G. Knepley 
1776858538eSMatthew G. Knepley   Level: developer
1786858538eSMatthew G. Knepley 
1796858538eSMatthew G. Knepley   Note: 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
1806858538eSMatthew G. Knepley 
1816858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMLocalizeCoordinate()`
1826858538eSMatthew G. Knepley */
1836858538eSMatthew G. Knepley PetscErrorCode DMLocalizeAddCoordinate_Internal(DM dm, PetscInt dim, const PetscScalar anchor[], const PetscScalar in[], PetscScalar out[])
1846858538eSMatthew G. Knepley {
1856858538eSMatthew G. Knepley   PetscInt d;
1866858538eSMatthew G. Knepley 
1876858538eSMatthew G. Knepley   PetscFunctionBegin;
1886858538eSMatthew G. Knepley   if (!dm->maxCell) {
1896858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) out[d] += in[d];
1906858538eSMatthew G. Knepley   } else {
1916858538eSMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1926858538eSMatthew G. Knepley       const PetscReal maxC = dm->maxCell[d];
1936858538eSMatthew G. Knepley 
1946858538eSMatthew G. Knepley       if ((dm->L[d] > 0.0) && (PetscAbsScalar(anchor[d] - in[d]) > maxC)) {
1956858538eSMatthew G. Knepley         const PetscScalar newCoord = PetscRealPart(anchor[d]) > PetscRealPart(in[d]) ? dm->L[d] + in[d] : in[d] - dm->L[d];
1966858538eSMatthew G. Knepley 
1976858538eSMatthew G. Knepley         if (PetscAbsScalar(newCoord - anchor[d]) > maxC)
1986858538eSMatthew 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]));
1996858538eSMatthew G. Knepley         out[d] += newCoord;
2006858538eSMatthew G. Knepley       } else {
2016858538eSMatthew G. Knepley         out[d] += in[d];
2026858538eSMatthew G. Knepley       }
2036858538eSMatthew G. Knepley     }
2046858538eSMatthew G. Knepley   }
2056858538eSMatthew G. Knepley   PetscFunctionReturn(0);
2066858538eSMatthew G. Knepley }
2076858538eSMatthew G. Knepley 
2086858538eSMatthew G. Knepley /*@
2096858538eSMatthew G. Knepley   DMGetCoordinatesLocalizedLocal - Check if the DM coordinates have been localized for cells on this process
2106858538eSMatthew G. Knepley 
2116858538eSMatthew G. Knepley   Not collective
2126858538eSMatthew G. Knepley 
2136858538eSMatthew G. Knepley   Input Parameter:
2146858538eSMatthew G. Knepley . dm - The DM
2156858538eSMatthew G. Knepley 
2166858538eSMatthew G. Knepley   Output Parameter:
2176858538eSMatthew G. Knepley   areLocalized - True if localized
2186858538eSMatthew G. Knepley 
2196858538eSMatthew G. Knepley   Level: developer
2206858538eSMatthew G. Knepley 
2216858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMGetCoordinatesLocalized()`, `DMSetPeriodicity()`
2226858538eSMatthew G. Knepley @*/
2236858538eSMatthew G. Knepley PetscErrorCode DMGetCoordinatesLocalizedLocal(DM dm, PetscBool *areLocalized)
2246858538eSMatthew G. Knepley {
2256858538eSMatthew G. Knepley   PetscFunctionBegin;
2266858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2276858538eSMatthew G. Knepley   PetscValidBoolPointer(areLocalized, 2);
2286858538eSMatthew G. Knepley   *areLocalized = dm->coordinates[1].dim < 0 ? PETSC_FALSE : PETSC_TRUE;
2296858538eSMatthew G. Knepley   PetscFunctionReturn(0);
2306858538eSMatthew G. Knepley }
2316858538eSMatthew G. Knepley 
2326858538eSMatthew G. Knepley /*@
2336858538eSMatthew G. Knepley   DMGetCoordinatesLocalized - Check if the DM coordinates have been localized for cells
2346858538eSMatthew G. Knepley 
2356858538eSMatthew G. Knepley   Collective on dm
2366858538eSMatthew G. Knepley 
2376858538eSMatthew G. Knepley   Input Parameter:
2386858538eSMatthew G. Knepley . dm - The DM
2396858538eSMatthew G. Knepley 
2406858538eSMatthew G. Knepley   Output Parameter:
2416858538eSMatthew G. Knepley   areLocalized - True if localized
2426858538eSMatthew G. Knepley 
2436858538eSMatthew G. Knepley   Level: developer
2446858538eSMatthew G. Knepley 
2456858538eSMatthew G. Knepley .seealso: `DMLocalizeCoordinates()`, `DMSetPeriodicity()`, `DMGetCoordinatesLocalizedLocal()`
2466858538eSMatthew G. Knepley @*/
2476858538eSMatthew G. Knepley PetscErrorCode DMGetCoordinatesLocalized(DM dm,PetscBool *areLocalized)
2486858538eSMatthew G. Knepley {
2496858538eSMatthew G. Knepley   PetscBool localized;
2506858538eSMatthew G. Knepley 
2516858538eSMatthew G. Knepley   PetscFunctionBegin;
2526858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2536858538eSMatthew G. Knepley   PetscValidBoolPointer(areLocalized, 2);
2546858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalizedLocal(dm, &localized));
2556858538eSMatthew G. Knepley   PetscCall(MPIU_Allreduce(&localized, areLocalized, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject) dm)));
2566858538eSMatthew G. Knepley   PetscFunctionReturn(0);
2576858538eSMatthew G. Knepley }
2586858538eSMatthew G. Knepley 
2596858538eSMatthew G. Knepley /*@
2606858538eSMatthew G. Knepley   DMLocalizeCoordinates - If a mesh is periodic, create local coordinates for cells having periodic faces
2616858538eSMatthew G. Knepley 
2626858538eSMatthew G. Knepley   Collective on dm
2636858538eSMatthew G. Knepley 
2646858538eSMatthew G. Knepley   Input Parameter:
2656858538eSMatthew G. Knepley . dm - The DM
2666858538eSMatthew G. Knepley 
2676858538eSMatthew G. Knepley   Level: developer
2686858538eSMatthew G. Knepley 
2696858538eSMatthew G. Knepley .seealso: `DMSetPeriodicity()`, `DMLocalizeCoordinate()`, `DMLocalizeAddCoordinate()`
2706858538eSMatthew G. Knepley @*/
2716858538eSMatthew G. Knepley PetscErrorCode DMLocalizeCoordinates(DM dm)
2726858538eSMatthew G. Knepley {
2736858538eSMatthew G. Knepley   DM               cdm, cdgdm, cplex, plex;
2746858538eSMatthew G. Knepley   PetscSection     cs, csDG;
2756858538eSMatthew G. Knepley   Vec              coordinates, cVec;
2766858538eSMatthew G. Knepley   PetscScalar     *coordsDG, *anchor, *localized;
2776858538eSMatthew G. Knepley   const PetscReal *L;
2786858538eSMatthew G. Knepley   PetscInt         Nc, vStart, vEnd, sStart, sEnd, newStart = PETSC_MAX_INT, newEnd = PETSC_MIN_INT, bs, coordSize;
2796858538eSMatthew G. Knepley   PetscBool        isLocalized, sparseLocalize = dm->sparseLocalize, useDG = PETSC_FALSE, useDGGlobal;
2806858538eSMatthew G. Knepley   PetscInt         maxHeight = 0, h;
2816858538eSMatthew G. Knepley   PetscInt        *pStart = NULL, *pEnd = NULL;
2826858538eSMatthew G. Knepley   MPI_Comm         comm;
2836858538eSMatthew G. Knepley 
2846858538eSMatthew G. Knepley   PetscFunctionBegin;
2856858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2866858538eSMatthew G. Knepley   PetscCall(DMGetPeriodicity(dm, NULL, &L));
2876858538eSMatthew G. Knepley   /* Cannot automatically localize without L and maxCell right now */
2886858538eSMatthew G. Knepley   if (!L) PetscFunctionReturn(0);
2896858538eSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
2906858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalized(dm, &isLocalized));
2916858538eSMatthew G. Knepley   if (isLocalized) PetscFunctionReturn(0);
2926858538eSMatthew G. Knepley 
2936858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2946858538eSMatthew G. Knepley   PetscCall(DMConvert(dm,  DMPLEX, &plex));
2956858538eSMatthew G. Knepley   PetscCall(DMConvert(cdm, DMPLEX, &cplex));
2966858538eSMatthew G. Knepley   if (cplex) {
2976858538eSMatthew G. Knepley     PetscCall(DMPlexGetDepthStratum(cplex, 0, &vStart, &vEnd));
2986858538eSMatthew G. Knepley     PetscCall(DMPlexGetMaxProjectionHeight(cplex, &maxHeight));
2996858538eSMatthew G. Knepley     PetscCall(DMGetWorkArray(dm, 2*(maxHeight + 1), MPIU_INT, &pStart));
3006858538eSMatthew G. Knepley     pEnd     = &pStart[maxHeight + 1];
3016858538eSMatthew G. Knepley     newStart = vStart;
3026858538eSMatthew G. Knepley     newEnd   = vEnd;
3036858538eSMatthew G. Knepley     for (h = 0; h <= maxHeight; h++) {
3046858538eSMatthew G. Knepley       PetscCall(DMPlexGetHeightStratum(cplex, h, &pStart[h], &pEnd[h]));
3056858538eSMatthew G. Knepley       newStart = PetscMin(newStart, pStart[h]);
3066858538eSMatthew G. Knepley       newEnd   = PetscMax(newEnd, pEnd[h]);
3076858538eSMatthew G. Knepley     }
3086858538eSMatthew G. Knepley   } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Coordinate localization requires a DMPLEX coordinate DM");
3096858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
3106858538eSMatthew G. Knepley   PetscCheck(coordinates, comm, PETSC_ERR_SUP, "Missing local coordinate vector");
3116858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateSection(dm, &cs));
3126858538eSMatthew G. Knepley   PetscCall(VecGetBlockSize(coordinates, &bs));
3136858538eSMatthew G. Knepley   PetscCall(PetscSectionGetChart(cs, &sStart, &sEnd));
3146858538eSMatthew G. Knepley 
3156858538eSMatthew G. Knepley   PetscCall(PetscSectionCreate(comm, &csDG));
3166858538eSMatthew G. Knepley   PetscCall(PetscSectionSetNumFields(csDG, 1));
3176858538eSMatthew G. Knepley   PetscCall(PetscSectionGetFieldComponents(cs, 0, &Nc));
3186858538eSMatthew G. Knepley   PetscCall(PetscSectionSetFieldComponents(csDG, 0, Nc));
3196858538eSMatthew G. Knepley   PetscCall(PetscSectionSetChart(csDG, newStart, newEnd));
3206858538eSMatthew G. Knepley   PetscCheck(bs == Nc, comm, PETSC_ERR_ARG_INCOMP, "Coordinate block size %" PetscInt_FMT " != %" PetscInt_FMT " number of components", bs, Nc);
3216858538eSMatthew G. Knepley 
3226858538eSMatthew G. Knepley   PetscCall(DMGetWorkArray(dm, 2 * Nc, MPIU_SCALAR, &anchor));
3236858538eSMatthew G. Knepley   localized = &anchor[Nc];
3246858538eSMatthew G. Knepley   for (h = 0; h <= maxHeight; h++) {
3256858538eSMatthew G. Knepley     PetscInt cStart = pStart[h], cEnd = pEnd[h], c;
3266858538eSMatthew G. Knepley 
3276858538eSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
3286858538eSMatthew G. Knepley       PetscScalar   *cellCoords = NULL;
3296858538eSMatthew G. Knepley       DMPolytopeType ct;
3306858538eSMatthew G. Knepley       PetscInt       dof, d, p;
3316858538eSMatthew G. Knepley 
3326858538eSMatthew G. Knepley       PetscCall(DMPlexGetCellType(plex, c, &ct));
3336858538eSMatthew G. Knepley       if (ct == DM_POLYTOPE_FV_GHOST) continue;
3346858538eSMatthew G. Knepley       PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
3356858538eSMatthew 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);
3366858538eSMatthew G. Knepley       for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[d];
3376858538eSMatthew G. Knepley       for (p = 0; p < dof/Nc; ++p) {
3386858538eSMatthew G. Knepley         PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p*Nc], localized));
3396858538eSMatthew G. Knepley         for (d = 0; d < Nc; ++d) if (cellCoords[p*Nc + d] != localized[d]) break;
3406858538eSMatthew G. Knepley         if (d < Nc) break;
3416858538eSMatthew G. Knepley       }
3426858538eSMatthew G. Knepley       if (p < dof/Nc) useDG = PETSC_TRUE;
3436858538eSMatthew G. Knepley       if (p < dof/Nc || !sparseLocalize) {
3446858538eSMatthew G. Knepley         PetscCall(PetscSectionSetDof(csDG, c, dof));
3456858538eSMatthew G. Knepley         PetscCall(PetscSectionSetFieldDof(csDG, c, 0, dof));
3466858538eSMatthew G. Knepley       }
3476858538eSMatthew G. Knepley       PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
3486858538eSMatthew G. Knepley     }
3496858538eSMatthew G. Knepley   }
3506858538eSMatthew G. Knepley   PetscCallMPI(MPI_Allreduce(&useDG, &useDGGlobal, 1, MPIU_BOOL, MPI_LOR, comm));
3516858538eSMatthew G. Knepley   if (!useDGGlobal) goto end;
3526858538eSMatthew G. Knepley 
3536858538eSMatthew G. Knepley   PetscCall(PetscSectionSetUp(csDG));
3546858538eSMatthew G. Knepley   PetscCall(PetscSectionGetStorageSize(csDG, &coordSize));
3556858538eSMatthew G. Knepley   PetscCall(VecCreate(PETSC_COMM_SELF, &cVec));
3566858538eSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject) cVec, "coordinates"));
3576858538eSMatthew G. Knepley   PetscCall(VecSetBlockSize(cVec, bs));
3586858538eSMatthew G. Knepley   PetscCall(VecSetSizes(cVec, coordSize, PETSC_DETERMINE));
3596858538eSMatthew G. Knepley   PetscCall(VecSetType(cVec, VECSTANDARD));
3606858538eSMatthew G. Knepley   PetscCall(VecGetArray(cVec, &coordsDG));
3616858538eSMatthew G. Knepley   for (h = 0; h <= maxHeight; h++) {
3626858538eSMatthew G. Knepley     PetscInt cStart = pStart[h], cEnd = pEnd[h], c;
3636858538eSMatthew G. Knepley 
3646858538eSMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
3656858538eSMatthew G. Knepley       PetscScalar *cellCoords = NULL;
3666858538eSMatthew G. Knepley       PetscInt     p = 0, q, dof, cdof, d, offDG;
3676858538eSMatthew G. Knepley 
3686858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(csDG, c, &cdof));
3696858538eSMatthew G. Knepley       if (!cdof) continue;
3706858538eSMatthew G. Knepley       PetscCall(DMPlexVecGetClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
3716858538eSMatthew G. Knepley       PetscCall(PetscSectionGetOffset(csDG, c, &offDG));
3726858538eSMatthew G. Knepley       // We need the cell to fit into [0, L]
373*7c48043bSMatthew G. Knepley       //   TODO The coordinates are set in closure order, which might not be the tensor order
3746858538eSMatthew G. Knepley       for (q = 0; q < dof/Nc; ++q) {
3756858538eSMatthew G. Knepley         // Select a trial anchor
3766858538eSMatthew G. Knepley         for (d = 0; d < Nc; ++d) anchor[d] = cellCoords[q*Nc+d];
3776858538eSMatthew G. Knepley         for (p = 0; p < dof/Nc; ++p) {
3786858538eSMatthew G. Knepley           PetscCall(DMLocalizeCoordinate_Internal(dm, Nc, anchor, &cellCoords[p*Nc], &coordsDG[offDG + p*Nc]));
3796858538eSMatthew G. Knepley           for (d = 0; d < Nc; ++d)
3806858538eSMatthew G. Knepley             if (L[d] > 0. && ((PetscRealPart(coordsDG[offDG + p*Nc + d]) < 0.) || (PetscRealPart(coordsDG[offDG + p*Nc + d]) > L[d]))) break;
3816858538eSMatthew G. Knepley           if (d < Nc) break;
3826858538eSMatthew G. Knepley         }
3836858538eSMatthew G. Knepley         if (p == dof/Nc) break;
3846858538eSMatthew G. Knepley       }
3856858538eSMatthew G. Knepley       PetscCheck(p == dof/Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " does not fit into the torus [0, L]", c);
3866858538eSMatthew G. Knepley       PetscCall(DMPlexVecRestoreClosure(cplex, cs, coordinates, c, &dof, &cellCoords));
3876858538eSMatthew G. Knepley     }
3886858538eSMatthew G. Knepley   }
3896858538eSMatthew G. Knepley   PetscCall(VecRestoreArray(cVec, &coordsDG));
3906858538eSMatthew G. Knepley   PetscCall(DMClone(cdm, &cdgdm));
3916858538eSMatthew G. Knepley   PetscCall(DMSetCellCoordinateDM(dm, cdgdm));
3926858538eSMatthew G. Knepley   PetscCall(DMSetCellCoordinateSection(dm, PETSC_DETERMINE, csDG));
3936858538eSMatthew G. Knepley   PetscCall(DMSetCellCoordinatesLocal(dm, cVec));
3946858538eSMatthew G. Knepley   PetscCall(VecDestroy(&cVec));
395*7c48043bSMatthew G. Knepley   // Convert the discretization
396*7c48043bSMatthew G. Knepley   {
397*7c48043bSMatthew G. Knepley     PetscFE         fe, dgfe;
398*7c48043bSMatthew G. Knepley     PetscSpace      P;
399*7c48043bSMatthew G. Knepley     PetscDualSpace  Q, dgQ;
400*7c48043bSMatthew G. Knepley     PetscQuadrature q, fq;
401*7c48043bSMatthew G. Knepley     PetscClassId    id;
402*7c48043bSMatthew G. Knepley 
403*7c48043bSMatthew G. Knepley     PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *) &fe));
404*7c48043bSMatthew G. Knepley     PetscCall(PetscObjectGetClassId((PetscObject) fe, &id));
405*7c48043bSMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
406*7c48043bSMatthew G. Knepley       PetscCall(PetscFEGetBasisSpace(fe, &P));
407*7c48043bSMatthew G. Knepley       PetscCall(PetscObjectReference((PetscObject) P));
408*7c48043bSMatthew G. Knepley       PetscCall(PetscFEGetDualSpace(fe, &Q));
409*7c48043bSMatthew G. Knepley       PetscCall(PetscDualSpaceDuplicate(Q, &dgQ));
410*7c48043bSMatthew G. Knepley       PetscCall(PetscDualSpaceLagrangeSetContinuity(dgQ, PETSC_FALSE));
411*7c48043bSMatthew G. Knepley       PetscCall(PetscDualSpaceSetUp(dgQ));
412*7c48043bSMatthew G. Knepley       PetscCall(PetscFEGetQuadrature(fe, &q));
413*7c48043bSMatthew G. Knepley       PetscCall(PetscObjectReference((PetscObject) q));
414*7c48043bSMatthew G. Knepley       PetscCall(PetscFEGetFaceQuadrature(fe, &fq));
415*7c48043bSMatthew G. Knepley       PetscCall(PetscObjectReference((PetscObject) fq));
416*7c48043bSMatthew G. Knepley       PetscCall(PetscFECreateFromSpaces(P, dgQ, q, fq, &dgfe));
417*7c48043bSMatthew G. Knepley       PetscCall(DMSetField(cdgdm, 0, NULL, (PetscObject) dgfe));
418*7c48043bSMatthew G. Knepley       PetscCall(PetscFEDestroy(&dgfe));
419*7c48043bSMatthew G. Knepley       PetscCall(DMCreateDS(cdgdm));
420*7c48043bSMatthew G. Knepley     }
421*7c48043bSMatthew G. Knepley   }
4226858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cdgdm));
4236858538eSMatthew G. Knepley 
4246858538eSMatthew G. Knepley end:
4256858538eSMatthew G. Knepley   PetscCall(DMRestoreWorkArray(dm, 2 * bs, MPIU_SCALAR, &anchor));
4266858538eSMatthew G. Knepley   PetscCall(DMRestoreWorkArray(dm, 2*(maxHeight + 1), MPIU_INT, &pStart));
4276858538eSMatthew G. Knepley   PetscCall(PetscSectionDestroy(&csDG));
4286858538eSMatthew G. Knepley   PetscCall(DMDestroy(&plex));
4296858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cplex));
4306858538eSMatthew G. Knepley   PetscFunctionReturn(0);
4316858538eSMatthew G. Knepley }
432