xref: /petsc/src/dm/interface/dmcoordinates.c (revision 356ea6bc55263ef03692b27b17f9e29fe78f4084)
16858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
26858538eSMatthew G. Knepley 
3e44f6aebSMatthew G. Knepley #include <petscdmplex.h> /* For DMCreateAffineCoordinates_Internal() */
46858538eSMatthew G. Knepley #include <petscsf.h>     /* For DMLocatePoints() */
56858538eSMatthew G. Knepley 
6d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, void *ctx)
7d71ae5a4SJacob Faibussowitsch {
86858538eSMatthew G. Knepley   DM  dm_coord, dmc_coord;
96858538eSMatthew G. Knepley   Vec coords, ccoords;
106858538eSMatthew G. Knepley   Mat inject;
114d86920dSPierre Jolivet 
126858538eSMatthew G. Knepley   PetscFunctionBegin;
136858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
146858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dmc, &dmc_coord));
156858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
166858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dmc, &ccoords));
176858538eSMatthew G. Knepley   if (coords && !ccoords) {
186858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords));
196858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
206858538eSMatthew G. Knepley     PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject));
216858538eSMatthew G. Knepley     PetscCall(MatRestrict(inject, coords, ccoords));
226858538eSMatthew G. Knepley     PetscCall(MatDestroy(&inject));
236858538eSMatthew G. Knepley     PetscCall(DMSetCoordinates(dmc, ccoords));
246858538eSMatthew G. Knepley     PetscCall(VecDestroy(&ccoords));
256858538eSMatthew G. Knepley   }
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
276858538eSMatthew G. Knepley }
286858538eSMatthew G. Knepley 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx)
30d71ae5a4SJacob Faibussowitsch {
316858538eSMatthew G. Knepley   DM          dm_coord, subdm_coord;
326858538eSMatthew G. Knepley   Vec         coords, ccoords, clcoords;
336858538eSMatthew G. Knepley   VecScatter *scat_i, *scat_g;
344d86920dSPierre Jolivet 
356858538eSMatthew G. Knepley   PetscFunctionBegin;
366858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
376858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(subdm, &subdm_coord));
386858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
396858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(subdm, &ccoords));
406858538eSMatthew G. Knepley   if (coords && !ccoords) {
416858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords));
426858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
436858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(subdm_coord, &clcoords));
446858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates"));
456858538eSMatthew G. Knepley     PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g));
466858538eSMatthew G. Knepley     PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
476858538eSMatthew G. Knepley     PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
486858538eSMatthew G. Knepley     PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
496858538eSMatthew G. Knepley     PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
506858538eSMatthew G. Knepley     PetscCall(DMSetCoordinates(subdm, ccoords));
516858538eSMatthew G. Knepley     PetscCall(DMSetCoordinatesLocal(subdm, clcoords));
526858538eSMatthew G. Knepley     PetscCall(VecScatterDestroy(&scat_i[0]));
536858538eSMatthew G. Knepley     PetscCall(VecScatterDestroy(&scat_g[0]));
546858538eSMatthew G. Knepley     PetscCall(VecDestroy(&ccoords));
556858538eSMatthew G. Knepley     PetscCall(VecDestroy(&clcoords));
566858538eSMatthew G. Knepley     PetscCall(PetscFree(scat_i));
576858538eSMatthew G. Knepley     PetscCall(PetscFree(scat_g));
586858538eSMatthew G. Knepley   }
593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
606858538eSMatthew G. Knepley }
616858538eSMatthew G. Knepley 
626858538eSMatthew G. Knepley /*@
63dce8aebaSBarry Smith   DMGetCoordinateDM - Gets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
646858538eSMatthew G. Knepley 
6520f4b53cSBarry Smith   Collective
666858538eSMatthew G. Knepley 
676858538eSMatthew G. Knepley   Input Parameter:
68dce8aebaSBarry Smith . dm - the `DM`
696858538eSMatthew G. Knepley 
706858538eSMatthew G. Knepley   Output Parameter:
71dce8aebaSBarry Smith . cdm - coordinate `DM`
726858538eSMatthew G. Knepley 
736858538eSMatthew G. Knepley   Level: intermediate
746858538eSMatthew G. Knepley 
75dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGSetCellCoordinateDM()`,
7660225df5SJacob Faibussowitsch 
776858538eSMatthew G. Knepley @*/
78d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
79d71ae5a4SJacob Faibussowitsch {
806858538eSMatthew G. Knepley   PetscFunctionBegin;
816858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
824f572ea9SToby Isaac   PetscAssertPointer(cdm, 2);
836858538eSMatthew G. Knepley   if (!dm->coordinates[0].dm) {
846858538eSMatthew G. Knepley     DM cdm;
856858538eSMatthew G. Knepley 
86dbbe0bcdSBarry Smith     PetscUseTypeMethod(dm, createcoordinatedm, &cdm);
876858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM"));
886858538eSMatthew G. Knepley     /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup
896858538eSMatthew G. Knepley      * until the call to CreateCoordinateDM) */
906858538eSMatthew G. Knepley     PetscCall(DMDestroy(&dm->coordinates[0].dm));
916858538eSMatthew G. Knepley     dm->coordinates[0].dm = cdm;
926858538eSMatthew G. Knepley   }
936858538eSMatthew G. Knepley   *cdm = dm->coordinates[0].dm;
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
956858538eSMatthew G. Knepley }
966858538eSMatthew G. Knepley 
976858538eSMatthew G. Knepley /*@
98dce8aebaSBarry Smith   DMSetCoordinateDM - Sets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
996858538eSMatthew G. Knepley 
10020f4b53cSBarry Smith   Logically Collective
1016858538eSMatthew G. Knepley 
1026858538eSMatthew G. Knepley   Input Parameters:
103dce8aebaSBarry Smith + dm  - the `DM`
104dce8aebaSBarry Smith - cdm - coordinate `DM`
1056858538eSMatthew G. Knepley 
1066858538eSMatthew G. Knepley   Level: intermediate
1076858538eSMatthew G. Knepley 
108dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMGetCellCoordinateDM()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`,
109dce8aebaSBarry Smith           `DMGSetCellCoordinateDM()`
1106858538eSMatthew G. Knepley @*/
111d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm)
112d71ae5a4SJacob Faibussowitsch {
1136858538eSMatthew G. Knepley   PetscFunctionBegin;
1146858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11590612307SJed Brown   if (cdm) PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
1166858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)cdm));
1176858538eSMatthew G. Knepley   PetscCall(DMDestroy(&dm->coordinates[0].dm));
1186858538eSMatthew G. Knepley   dm->coordinates[0].dm = cdm;
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1206858538eSMatthew G. Knepley }
1216858538eSMatthew G. Knepley 
1226858538eSMatthew G. Knepley /*@
123dce8aebaSBarry Smith   DMGetCellCoordinateDM - Gets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1246858538eSMatthew G. Knepley 
12520f4b53cSBarry Smith   Collective
1266858538eSMatthew G. Knepley 
1276858538eSMatthew G. Knepley   Input Parameter:
128dce8aebaSBarry Smith . dm - the `DM`
1296858538eSMatthew G. Knepley 
1306858538eSMatthew G. Knepley   Output Parameter:
1310b3275a6SBarry Smith . cdm - cellwise coordinate `DM`, or `NULL` if they are not defined
1326858538eSMatthew G. Knepley 
1336858538eSMatthew G. Knepley   Level: intermediate
1346858538eSMatthew G. Knepley 
135dce8aebaSBarry Smith   Note:
136dce8aebaSBarry Smith   Call `DMLocalizeCoordinates()` to automatically create cellwise coordinates for periodic geometries.
137dce8aebaSBarry Smith 
138dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
139dce8aebaSBarry Smith           `DMLocalizeCoordinates()`, `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
1406858538eSMatthew G. Knepley @*/
141d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm)
142d71ae5a4SJacob Faibussowitsch {
1436858538eSMatthew G. Knepley   PetscFunctionBegin;
1446858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1454f572ea9SToby Isaac   PetscAssertPointer(cdm, 2);
1466858538eSMatthew G. Knepley   *cdm = dm->coordinates[1].dm;
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1486858538eSMatthew G. Knepley }
1496858538eSMatthew G. Knepley 
1506858538eSMatthew G. Knepley /*@
151dce8aebaSBarry Smith   DMSetCellCoordinateDM - Sets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1526858538eSMatthew G. Knepley 
15320f4b53cSBarry Smith   Logically Collective
1546858538eSMatthew G. Knepley 
1556858538eSMatthew G. Knepley   Input Parameters:
156dce8aebaSBarry Smith + dm  - the `DM`
157dce8aebaSBarry Smith - cdm - cellwise coordinate `DM`
1586858538eSMatthew G. Knepley 
1596858538eSMatthew G. Knepley   Level: intermediate
1606858538eSMatthew G. Knepley 
161dce8aebaSBarry Smith   Note:
16235cb6cd3SPierre Jolivet   As opposed to `DMSetCoordinateDM()` these coordinates are useful for discontinuous Galerkin methods since they support coordinate fields that are discontinuous at cell boundaries.
163dce8aebaSBarry Smith 
164dce8aebaSBarry Smith .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
165dce8aebaSBarry Smith           `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
1666858538eSMatthew G. Knepley @*/
167d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm)
168d71ae5a4SJacob Faibussowitsch {
1696858538eSMatthew G. Knepley   PetscInt dim;
1706858538eSMatthew G. Knepley 
1716858538eSMatthew G. Knepley   PetscFunctionBegin;
1726858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1736858538eSMatthew G. Knepley   if (cdm) {
1746858538eSMatthew G. Knepley     PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
1756858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &dim));
1766858538eSMatthew G. Knepley     dm->coordinates[1].dim = dim;
1776858538eSMatthew G. Knepley   }
1786858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)cdm));
1796858538eSMatthew G. Knepley   PetscCall(DMDestroy(&dm->coordinates[1].dm));
1806858538eSMatthew G. Knepley   dm->coordinates[1].dm = cdm;
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1826858538eSMatthew G. Knepley }
1836858538eSMatthew G. Knepley 
1846858538eSMatthew G. Knepley /*@
1850b3275a6SBarry Smith   DMGetCoordinateDim - Retrieve the dimension of the embedding space for coordinate values. For example a mesh on the surface of a sphere would have a 3 dimensional embedding space
1866858538eSMatthew G. Knepley 
1876858538eSMatthew G. Knepley   Not Collective
1886858538eSMatthew G. Knepley 
1896858538eSMatthew G. Knepley   Input Parameter:
190dce8aebaSBarry Smith . dm - The `DM` object
1916858538eSMatthew G. Knepley 
1926858538eSMatthew G. Knepley   Output Parameter:
1936858538eSMatthew G. Knepley . dim - The embedding dimension
1946858538eSMatthew G. Knepley 
1956858538eSMatthew G. Knepley   Level: intermediate
1966858538eSMatthew G. Knepley 
197dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
1986858538eSMatthew G. Knepley @*/
199d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim)
200d71ae5a4SJacob Faibussowitsch {
2016858538eSMatthew G. Knepley   PetscFunctionBegin;
2026858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2034f572ea9SToby Isaac   PetscAssertPointer(dim, 2);
2046858538eSMatthew G. Knepley   if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim;
2056858538eSMatthew G. Knepley   *dim = dm->coordinates[0].dim;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2076858538eSMatthew G. Knepley }
2086858538eSMatthew G. Knepley 
2096858538eSMatthew G. Knepley /*@
2106858538eSMatthew G. Knepley   DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values.
2116858538eSMatthew G. Knepley 
2126858538eSMatthew G. Knepley   Not Collective
2136858538eSMatthew G. Knepley 
2146858538eSMatthew G. Knepley   Input Parameters:
215dce8aebaSBarry Smith + dm  - The `DM` object
2166858538eSMatthew G. Knepley - dim - The embedding dimension
2176858538eSMatthew G. Knepley 
2186858538eSMatthew G. Knepley   Level: intermediate
2196858538eSMatthew G. Knepley 
220dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2216858538eSMatthew G. Knepley @*/
222d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim)
223d71ae5a4SJacob Faibussowitsch {
2246858538eSMatthew G. Knepley   PetscDS  ds;
2256858538eSMatthew G. Knepley   PetscInt Nds, n;
2266858538eSMatthew G. Knepley 
2276858538eSMatthew G. Knepley   PetscFunctionBegin;
2286858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2296858538eSMatthew G. Knepley   dm->coordinates[0].dim = dim;
2306858538eSMatthew G. Knepley   if (dm->dim >= 0) {
2316858538eSMatthew G. Knepley     PetscCall(DMGetNumDS(dm, &Nds));
2326858538eSMatthew G. Knepley     for (n = 0; n < Nds; ++n) {
23307218a29SMatthew G. Knepley       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
2346858538eSMatthew G. Knepley       PetscCall(PetscDSSetCoordinateDimension(ds, dim));
2356858538eSMatthew G. Knepley     }
2366858538eSMatthew G. Knepley   }
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2386858538eSMatthew G. Knepley }
2396858538eSMatthew G. Knepley 
2406858538eSMatthew G. Knepley /*@
2410b3275a6SBarry Smith   DMGetCoordinateSection - Retrieve the `PetscSection` of coordinate values over the mesh.
2426858538eSMatthew G. Knepley 
24320f4b53cSBarry Smith   Collective
2446858538eSMatthew G. Knepley 
2456858538eSMatthew G. Knepley   Input Parameter:
246dce8aebaSBarry Smith . dm - The `DM` object
2476858538eSMatthew G. Knepley 
2486858538eSMatthew G. Knepley   Output Parameter:
249dce8aebaSBarry Smith . section - The `PetscSection` object
2506858538eSMatthew G. Knepley 
2516858538eSMatthew G. Knepley   Level: intermediate
2526858538eSMatthew G. Knepley 
253dce8aebaSBarry Smith   Note:
254dce8aebaSBarry Smith   This just retrieves the local section from the coordinate `DM`. In other words,
255dce8aebaSBarry Smith .vb
256dce8aebaSBarry Smith   DMGetCoordinateDM(dm, &cdm);
257dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
258dce8aebaSBarry Smith .ve
259dce8aebaSBarry Smith 
2600b3275a6SBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2616858538eSMatthew G. Knepley @*/
262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section)
263d71ae5a4SJacob Faibussowitsch {
2646858538eSMatthew G. Knepley   DM cdm;
2656858538eSMatthew G. Knepley 
2666858538eSMatthew G. Knepley   PetscFunctionBegin;
2676858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2684f572ea9SToby Isaac   PetscAssertPointer(section, 2);
2696858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2706858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, section));
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2726858538eSMatthew G. Knepley }
2736858538eSMatthew G. Knepley 
2746858538eSMatthew G. Knepley /*@
2750b3275a6SBarry Smith   DMSetCoordinateSection - Set the `PetscSection` of coordinate values over the mesh.
2766858538eSMatthew G. Knepley 
2776858538eSMatthew G. Knepley   Not Collective
2786858538eSMatthew G. Knepley 
2796858538eSMatthew G. Knepley   Input Parameters:
280dce8aebaSBarry Smith + dm      - The `DM` object
281dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
282dce8aebaSBarry Smith - section - The `PetscSection` object
2836858538eSMatthew G. Knepley 
2846858538eSMatthew G. Knepley   Level: intermediate
2856858538eSMatthew G. Knepley 
286dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2876858538eSMatthew G. Knepley @*/
288d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section)
289d71ae5a4SJacob Faibussowitsch {
2906858538eSMatthew G. Knepley   DM cdm;
2916858538eSMatthew G. Knepley 
2926858538eSMatthew G. Knepley   PetscFunctionBegin;
2936858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2946858538eSMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
2956858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2966858538eSMatthew G. Knepley   PetscCall(DMSetLocalSection(cdm, section));
2976858538eSMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
2986858538eSMatthew G. Knepley     PetscInt d = PETSC_DEFAULT;
2996858538eSMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
3006858538eSMatthew G. Knepley 
3016858538eSMatthew G. Knepley     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3026858538eSMatthew G. Knepley     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
3036858538eSMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
3046858538eSMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
3056858538eSMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
3066858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(section, v, &dd));
3079371c9d4SSatish Balay       if (dd) {
3089371c9d4SSatish Balay         d = dd;
3099371c9d4SSatish Balay         break;
3109371c9d4SSatish Balay       }
3116858538eSMatthew G. Knepley     }
3126858538eSMatthew G. Knepley     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
3136858538eSMatthew G. Knepley   }
3143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3156858538eSMatthew G. Knepley }
3166858538eSMatthew G. Knepley 
3176858538eSMatthew G. Knepley /*@
3180b3275a6SBarry Smith   DMGetCellCoordinateSection - Retrieve the `PetscSection` of cellwise coordinate values over the mesh.
3196858538eSMatthew G. Knepley 
32020f4b53cSBarry Smith   Collective
3216858538eSMatthew G. Knepley 
3226858538eSMatthew G. Knepley   Input Parameter:
323dce8aebaSBarry Smith . dm - The `DM` object
3246858538eSMatthew G. Knepley 
3256858538eSMatthew G. Knepley   Output Parameter:
32620f4b53cSBarry Smith . section - The `PetscSection` object, or `NULL` if no cellwise coordinates are defined
3276858538eSMatthew G. Knepley 
3286858538eSMatthew G. Knepley   Level: intermediate
3296858538eSMatthew G. Knepley 
330dce8aebaSBarry Smith   Note:
331dce8aebaSBarry Smith   This just retrieves the local section from the cell coordinate `DM`. In other words,
332dce8aebaSBarry Smith .vb
333dce8aebaSBarry Smith   DMGetCellCoordinateDM(dm, &cdm);
334dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
335dce8aebaSBarry Smith .ve
336dce8aebaSBarry Smith 
337dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
3386858538eSMatthew G. Knepley @*/
339d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section)
340d71ae5a4SJacob Faibussowitsch {
3416858538eSMatthew G. Knepley   DM cdm;
3426858538eSMatthew G. Knepley 
3436858538eSMatthew G. Knepley   PetscFunctionBegin;
3446858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3454f572ea9SToby Isaac   PetscAssertPointer(section, 2);
3466858538eSMatthew G. Knepley   *section = NULL;
3476858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3486858538eSMatthew G. Knepley   if (cdm) PetscCall(DMGetLocalSection(cdm, section));
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3506858538eSMatthew G. Knepley }
3516858538eSMatthew G. Knepley 
3526858538eSMatthew G. Knepley /*@
3530b3275a6SBarry Smith   DMSetCellCoordinateSection - Set the `PetscSection` of cellwise coordinate values over the mesh.
3546858538eSMatthew G. Knepley 
3556858538eSMatthew G. Knepley   Not Collective
3566858538eSMatthew G. Knepley 
3576858538eSMatthew G. Knepley   Input Parameters:
358dce8aebaSBarry Smith + dm      - The `DM` object
359dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
360dce8aebaSBarry Smith - section - The `PetscSection` object for a cellwise layout
3616858538eSMatthew G. Knepley 
3626858538eSMatthew G. Knepley   Level: intermediate
3636858538eSMatthew G. Knepley 
364dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
3656858538eSMatthew G. Knepley @*/
366d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section)
367d71ae5a4SJacob Faibussowitsch {
3686858538eSMatthew G. Knepley   DM cdm;
3696858538eSMatthew G. Knepley 
3706858538eSMatthew G. Knepley   PetscFunctionBegin;
3716858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3726858538eSMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
3736858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3746858538eSMatthew G. Knepley   PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates");
3756858538eSMatthew G. Knepley   PetscCall(DMSetLocalSection(cdm, section));
3766858538eSMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
3776858538eSMatthew G. Knepley     PetscInt d = PETSC_DEFAULT;
3786858538eSMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
3796858538eSMatthew G. Knepley 
3806858538eSMatthew G. Knepley     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3816858538eSMatthew G. Knepley     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
3826858538eSMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
3836858538eSMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
3846858538eSMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
3856858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(section, v, &dd));
3869371c9d4SSatish Balay       if (dd) {
3879371c9d4SSatish Balay         d = dd;
3889371c9d4SSatish Balay         break;
3899371c9d4SSatish Balay       }
3906858538eSMatthew G. Knepley     }
3916858538eSMatthew G. Knepley     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
3926858538eSMatthew G. Knepley   }
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3946858538eSMatthew G. Knepley }
3956858538eSMatthew G. Knepley 
3966858538eSMatthew G. Knepley /*@
397dce8aebaSBarry Smith   DMGetCoordinates - Gets a global vector with the coordinates associated with the `DM`.
3986858538eSMatthew G. Knepley 
399*356ea6bcSBarry Smith   Collective if the global vector with coordinates has not been set yet but the local vector with coordinates has been set
4006858538eSMatthew G. Knepley 
4016858538eSMatthew G. Knepley   Input Parameter:
402dce8aebaSBarry Smith . dm - the `DM`
4036858538eSMatthew G. Knepley 
4046858538eSMatthew G. Knepley   Output Parameter:
4056858538eSMatthew G. Knepley . c - global coordinate vector
4066858538eSMatthew G. Knepley 
407dce8aebaSBarry Smith   Level: intermediate
408dce8aebaSBarry Smith 
409dce8aebaSBarry Smith   Notes:
410dce8aebaSBarry Smith   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
4110b3275a6SBarry Smith   destroyed `c` will no longer be valid.
4126858538eSMatthew G. Knepley 
413*356ea6bcSBarry Smith   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates), see `DMGetCoordinatesLocal()`.
4146858538eSMatthew G. Knepley 
415dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
4166858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
4176858538eSMatthew G. Knepley 
418*356ea6bcSBarry Smith   Does not work for `DMSTAG`
419*356ea6bcSBarry Smith 
420dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
4216858538eSMatthew G. Knepley @*/
422d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
423d71ae5a4SJacob Faibussowitsch {
4246858538eSMatthew G. Knepley   PetscFunctionBegin;
4256858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4264f572ea9SToby Isaac   PetscAssertPointer(c, 2);
4276858538eSMatthew G. Knepley   if (!dm->coordinates[0].x && dm->coordinates[0].xl) {
4286858538eSMatthew G. Knepley     DM cdm = NULL;
4296858538eSMatthew G. Knepley 
4306858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
4316858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x));
4326858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates"));
4336858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
4346858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
4356858538eSMatthew G. Knepley   }
4366858538eSMatthew G. Knepley   *c = dm->coordinates[0].x;
4373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4386858538eSMatthew G. Knepley }
4396858538eSMatthew G. Knepley 
4406858538eSMatthew G. Knepley /*@
441dce8aebaSBarry Smith   DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates
4426858538eSMatthew G. Knepley 
443*356ea6bcSBarry Smith   Logically Collective
4446858538eSMatthew G. Knepley 
4456858538eSMatthew G. Knepley   Input Parameters:
446dce8aebaSBarry Smith + dm - the `DM`
4476858538eSMatthew G. Knepley - c  - coordinate vector
4486858538eSMatthew G. Knepley 
4496858538eSMatthew G. Knepley   Level: intermediate
4506858538eSMatthew G. Knepley 
451dce8aebaSBarry Smith   Notes:
452dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
453dce8aebaSBarry Smith 
45420f4b53cSBarry Smith   The vector `c` can be destroyed after the call
455dce8aebaSBarry Smith 
456dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
4576858538eSMatthew G. Knepley @*/
458d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinates(DM dm, Vec c)
459d71ae5a4SJacob Faibussowitsch {
4606858538eSMatthew G. Knepley   PetscFunctionBegin;
4616858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4626858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
4636858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
4646858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
4656858538eSMatthew G. Knepley   dm->coordinates[0].x = c;
4666858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
4676858538eSMatthew G. Knepley   PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL));
4686858538eSMatthew G. Knepley   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL));
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4706858538eSMatthew G. Knepley }
4716858538eSMatthew G. Knepley 
4726858538eSMatthew G. Knepley /*@
473dce8aebaSBarry Smith   DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`.
4746858538eSMatthew G. Knepley 
47520f4b53cSBarry Smith   Collective
4766858538eSMatthew G. Knepley 
4776858538eSMatthew G. Knepley   Input Parameter:
478dce8aebaSBarry Smith . dm - the `DM`
4796858538eSMatthew G. Knepley 
4806858538eSMatthew G. Knepley   Output Parameter:
4816858538eSMatthew G. Knepley . c - global coordinate vector
4826858538eSMatthew G. Knepley 
483dce8aebaSBarry Smith   Level: intermediate
484dce8aebaSBarry Smith 
485dce8aebaSBarry Smith   Notes:
486dce8aebaSBarry Smith   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
4870b3275a6SBarry Smith   destroyed `c` will no longer be valid.
4886858538eSMatthew G. Knepley 
4896858538eSMatthew G. Knepley   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
4906858538eSMatthew G. Knepley 
491dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
4926858538eSMatthew G. Knepley @*/
493d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c)
494d71ae5a4SJacob Faibussowitsch {
4956858538eSMatthew G. Knepley   PetscFunctionBegin;
4966858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4974f572ea9SToby Isaac   PetscAssertPointer(c, 2);
4986858538eSMatthew G. Knepley   if (!dm->coordinates[1].x && dm->coordinates[1].xl) {
4996858538eSMatthew G. Knepley     DM cdm = NULL;
5006858538eSMatthew G. Knepley 
50111ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
5026858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x));
5036858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates"));
5046858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
5056858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
5066858538eSMatthew G. Knepley   }
5076858538eSMatthew G. Knepley   *c = dm->coordinates[1].x;
5083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5096858538eSMatthew G. Knepley }
5106858538eSMatthew G. Knepley 
5116858538eSMatthew G. Knepley /*@
512dce8aebaSBarry Smith   DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates
5136858538eSMatthew G. Knepley 
51420f4b53cSBarry Smith   Collective
5156858538eSMatthew G. Knepley 
5166858538eSMatthew G. Knepley   Input Parameters:
517dce8aebaSBarry Smith + dm - the `DM`
5186858538eSMatthew G. Knepley - c  - cellwise coordinate vector
5196858538eSMatthew G. Knepley 
5206858538eSMatthew G. Knepley   Level: intermediate
5216858538eSMatthew G. Knepley 
522dce8aebaSBarry Smith   Notes:
523dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
524dce8aebaSBarry Smith 
52520f4b53cSBarry Smith   The vector `c` should be destroyed by the caller.
526dce8aebaSBarry Smith 
527dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
5286858538eSMatthew G. Knepley @*/
529d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinates(DM dm, Vec c)
530d71ae5a4SJacob Faibussowitsch {
5316858538eSMatthew G. Knepley   PetscFunctionBegin;
5326858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5336858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
5346858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
5356858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
5366858538eSMatthew G. Knepley   dm->coordinates[1].x = c;
5376858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
5383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5396858538eSMatthew G. Knepley }
5406858538eSMatthew G. Knepley 
5416858538eSMatthew G. Knepley /*@
542dce8aebaSBarry Smith   DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
5436858538eSMatthew G. Knepley 
54420f4b53cSBarry Smith   Collective
5456858538eSMatthew G. Knepley 
5466858538eSMatthew G. Knepley   Input Parameter:
547dce8aebaSBarry Smith . dm - the `DM`
5486858538eSMatthew G. Knepley 
5496858538eSMatthew G. Knepley   Level: advanced
5506858538eSMatthew G. Knepley 
551dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMGetCoordinatesLocalNoncollective()`
5526858538eSMatthew G. Knepley @*/
553d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm)
554d71ae5a4SJacob Faibussowitsch {
5556858538eSMatthew G. Knepley   PetscFunctionBegin;
5566858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5576858538eSMatthew G. Knepley   if (!dm->coordinates[0].xl && dm->coordinates[0].x) {
5586858538eSMatthew G. Knepley     DM       cdm = NULL;
5599d92f5ffSMatthew G. Knepley     PetscInt bs;
5606858538eSMatthew G. Knepley 
5616858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
5626858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
56321c42226SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "Local Coordinates"));
5649d92f5ffSMatthew G. Knepley     // If the size of the vector is 0, it will not get the right block size
5659d92f5ffSMatthew G. Knepley     PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs));
5669d92f5ffSMatthew G. Knepley     PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs));
5676858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
5686858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5696858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5706858538eSMatthew G. Knepley   }
5713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5726858538eSMatthew G. Knepley }
5736858538eSMatthew G. Knepley 
5746858538eSMatthew G. Knepley /*@
575dce8aebaSBarry Smith   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.
5766858538eSMatthew G. Knepley 
57720f4b53cSBarry Smith   Collective the first time it is called
5786858538eSMatthew G. Knepley 
5796858538eSMatthew G. Knepley   Input Parameter:
580dce8aebaSBarry Smith . dm - the `DM`
5816858538eSMatthew G. Knepley 
5826858538eSMatthew G. Knepley   Output Parameter:
5836858538eSMatthew G. Knepley . c - coordinate vector
5846858538eSMatthew G. Knepley 
585dce8aebaSBarry Smith   Level: intermediate
586dce8aebaSBarry Smith 
587dce8aebaSBarry Smith   Notes:
5880b3275a6SBarry Smith   This is a borrowed reference, so the user should NOT destroy `c`
5896858538eSMatthew G. Knepley 
5906858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
5916858538eSMatthew G. Knepley 
592dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
5936858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
5946858538eSMatthew G. Knepley 
595dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
5966858538eSMatthew G. Knepley @*/
597d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
598d71ae5a4SJacob Faibussowitsch {
5996858538eSMatthew G. Knepley   PetscFunctionBegin;
6006858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6014f572ea9SToby Isaac   PetscAssertPointer(c, 2);
6026858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
6036858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
6043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6056858538eSMatthew G. Knepley }
6066858538eSMatthew G. Knepley 
6076858538eSMatthew G. Knepley /*@
608dce8aebaSBarry Smith   DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.
6096858538eSMatthew G. Knepley 
61020f4b53cSBarry Smith   Not Collective
6116858538eSMatthew G. Knepley 
6126858538eSMatthew G. Knepley   Input Parameter:
613dce8aebaSBarry Smith . dm - the `DM`
6146858538eSMatthew G. Knepley 
6156858538eSMatthew G. Knepley   Output Parameter:
6166858538eSMatthew G. Knepley . c - coordinate vector
6176858538eSMatthew G. Knepley 
6186858538eSMatthew G. Knepley   Level: advanced
6196858538eSMatthew G. Knepley 
620dce8aebaSBarry Smith   Note:
621dce8aebaSBarry Smith   A previous call to  `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.
622dce8aebaSBarry Smith 
623dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6246858538eSMatthew G. Knepley @*/
625d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
626d71ae5a4SJacob Faibussowitsch {
6276858538eSMatthew G. Knepley   PetscFunctionBegin;
6286858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6294f572ea9SToby Isaac   PetscAssertPointer(c, 2);
6306858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
6316858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
6323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6336858538eSMatthew G. Knepley }
6346858538eSMatthew G. Knepley 
6356858538eSMatthew G. Knepley /*@
636dce8aebaSBarry Smith   DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.
6376858538eSMatthew G. Knepley 
63820f4b53cSBarry Smith   Not Collective
6396858538eSMatthew G. Knepley 
6406858538eSMatthew G. Knepley   Input Parameters:
641dce8aebaSBarry Smith + dm - the `DM`
642dce8aebaSBarry Smith - p  - the `IS` of points whose coordinates will be returned
6436858538eSMatthew G. Knepley 
6446858538eSMatthew G. Knepley   Output Parameters:
6450b3275a6SBarry Smith + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in `p`, and DOFs correspond to coordinates
6460b3275a6SBarry Smith - pCoord        - the `Vec` with coordinates of points in `p`
6476858538eSMatthew G. Knepley 
648dce8aebaSBarry Smith   Level: advanced
649dce8aebaSBarry Smith 
650dce8aebaSBarry Smith   Notes:
651dce8aebaSBarry Smith   `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.
6526858538eSMatthew G. Knepley 
6536858538eSMatthew G. Knepley   This creates a new vector, so the user SHOULD destroy this vector
6546858538eSMatthew G. Knepley 
6556858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
6566858538eSMatthew G. Knepley 
657dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
6586858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
6596858538eSMatthew G. Knepley 
660dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6616858538eSMatthew G. Knepley @*/
662d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
663d71ae5a4SJacob Faibussowitsch {
6646858538eSMatthew G. Knepley   DM                 cdm;
6656858538eSMatthew G. Knepley   PetscSection       cs, newcs;
6666858538eSMatthew G. Knepley   Vec                coords;
6676858538eSMatthew G. Knepley   const PetscScalar *arr;
6686858538eSMatthew G. Knepley   PetscScalar       *newarr = NULL;
6696858538eSMatthew G. Knepley   PetscInt           n;
6706858538eSMatthew G. Knepley 
6716858538eSMatthew G. Knepley   PetscFunctionBegin;
6726858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6736858538eSMatthew G. Knepley   PetscValidHeaderSpecific(p, IS_CLASSID, 2);
6744f572ea9SToby Isaac   if (pCoordSection) PetscAssertPointer(pCoordSection, 3);
6754f572ea9SToby Isaac   if (pCoord) PetscAssertPointer(pCoord, 4);
6766858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
6776858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, &cs));
6786858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
6796858538eSMatthew G. Knepley   PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
6806858538eSMatthew G. Knepley   PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
6816858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(coords, &arr));
6826858538eSMatthew G. Knepley   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
6836858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(coords, &arr));
6846858538eSMatthew G. Knepley   if (pCoord) {
6856858538eSMatthew G. Knepley     PetscCall(PetscSectionGetStorageSize(newcs, &n));
6866858538eSMatthew G. Knepley     /* set array in two steps to mimic PETSC_OWN_POINTER */
6876858538eSMatthew G. Knepley     PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
6886858538eSMatthew G. Knepley     PetscCall(VecReplaceArray(*pCoord, newarr));
6896858538eSMatthew G. Knepley   } else {
6906858538eSMatthew G. Knepley     PetscCall(PetscFree(newarr));
6916858538eSMatthew G. Knepley   }
6929371c9d4SSatish Balay   if (pCoordSection) {
6939371c9d4SSatish Balay     *pCoordSection = newcs;
6949371c9d4SSatish Balay   } else PetscCall(PetscSectionDestroy(&newcs));
6953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6966858538eSMatthew G. Knepley }
6976858538eSMatthew G. Knepley 
6986858538eSMatthew G. Knepley /*@
699dce8aebaSBarry Smith   DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates
7006858538eSMatthew G. Knepley 
70120f4b53cSBarry Smith   Not Collective
7026858538eSMatthew G. Knepley 
7036858538eSMatthew G. Knepley   Input Parameters:
704dce8aebaSBarry Smith + dm - the `DM`
7056858538eSMatthew G. Knepley - c  - coordinate vector
7066858538eSMatthew G. Knepley 
707dce8aebaSBarry Smith   Level: intermediate
708dce8aebaSBarry Smith 
7096858538eSMatthew G. Knepley   Notes:
710dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
711dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
7126858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
7136858538eSMatthew G. Knepley 
7140b3275a6SBarry Smith   The vector `c` should be destroyed by the caller.
7156858538eSMatthew G. Knepley 
716dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
7176858538eSMatthew G. Knepley @*/
718d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
719d71ae5a4SJacob Faibussowitsch {
7206858538eSMatthew G. Knepley   PetscFunctionBegin;
7216858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7226858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
7236858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
7246858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
7256858538eSMatthew G. Knepley   dm->coordinates[0].xl = c;
7266858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
7273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7286858538eSMatthew G. Knepley }
7296858538eSMatthew G. Knepley 
7306858538eSMatthew G. Knepley /*@
731dce8aebaSBarry Smith   DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
7326858538eSMatthew G. Knepley 
73320f4b53cSBarry Smith   Collective
7346858538eSMatthew G. Knepley 
7356858538eSMatthew G. Knepley   Input Parameter:
736dce8aebaSBarry Smith . dm - the `DM`
7376858538eSMatthew G. Knepley 
7386858538eSMatthew G. Knepley   Level: advanced
7396858538eSMatthew G. Knepley 
740dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
7416858538eSMatthew G. Knepley @*/
742d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
743d71ae5a4SJacob Faibussowitsch {
7446858538eSMatthew G. Knepley   PetscFunctionBegin;
7456858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7466858538eSMatthew G. Knepley   if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
7476858538eSMatthew G. Knepley     DM cdm = NULL;
7486858538eSMatthew G. Knepley 
74911ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
7506858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
7516858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
7526858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7536858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7546858538eSMatthew G. Knepley   }
7553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7566858538eSMatthew G. Knepley }
7576858538eSMatthew G. Knepley 
7586858538eSMatthew G. Knepley /*@
759dce8aebaSBarry Smith   DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.
7606858538eSMatthew G. Knepley 
76120f4b53cSBarry Smith   Collective
7626858538eSMatthew G. Knepley 
7636858538eSMatthew G. Knepley   Input Parameter:
764dce8aebaSBarry Smith . dm - the `DM`
7656858538eSMatthew G. Knepley 
7666858538eSMatthew G. Knepley   Output Parameter:
7676858538eSMatthew G. Knepley . c - coordinate vector
7686858538eSMatthew G. Knepley 
769dce8aebaSBarry Smith   Level: intermediate
770dce8aebaSBarry Smith 
771dce8aebaSBarry Smith   Notes:
7726858538eSMatthew G. Knepley   This is a borrowed reference, so the user should NOT destroy this vector
7736858538eSMatthew G. Knepley 
7746858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
7756858538eSMatthew G. Knepley 
776dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
7776858538eSMatthew G. Knepley @*/
778d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
779d71ae5a4SJacob Faibussowitsch {
7806858538eSMatthew G. Knepley   PetscFunctionBegin;
7816858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7824f572ea9SToby Isaac   PetscAssertPointer(c, 2);
7836858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
7846858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
7853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7866858538eSMatthew G. Knepley }
7876858538eSMatthew G. Knepley 
7886858538eSMatthew G. Knepley /*@
789dce8aebaSBarry Smith   DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.
7906858538eSMatthew G. Knepley 
79120f4b53cSBarry Smith   Not Collective
7926858538eSMatthew G. Knepley 
7936858538eSMatthew G. Knepley   Input Parameter:
794dce8aebaSBarry Smith . dm - the `DM`
7956858538eSMatthew G. Knepley 
7966858538eSMatthew G. Knepley   Output Parameter:
7976858538eSMatthew G. Knepley . c - cellwise coordinate vector
7986858538eSMatthew G. Knepley 
7996858538eSMatthew G. Knepley   Level: advanced
8006858538eSMatthew G. Knepley 
801dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
8026858538eSMatthew G. Knepley @*/
803d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
804d71ae5a4SJacob Faibussowitsch {
8056858538eSMatthew G. Knepley   PetscFunctionBegin;
8066858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8074f572ea9SToby Isaac   PetscAssertPointer(c, 2);
8086858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
8096858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
8103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8116858538eSMatthew G. Knepley }
8126858538eSMatthew G. Knepley 
8136858538eSMatthew G. Knepley /*@
814dce8aebaSBarry Smith   DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates
8156858538eSMatthew G. Knepley 
81620f4b53cSBarry Smith   Not Collective
8176858538eSMatthew G. Knepley 
8186858538eSMatthew G. Knepley   Input Parameters:
819dce8aebaSBarry Smith + dm - the `DM`
8206858538eSMatthew G. Knepley - c  - cellwise coordinate vector
8216858538eSMatthew G. Knepley 
822dce8aebaSBarry Smith   Level: intermediate
823dce8aebaSBarry Smith 
8246858538eSMatthew G. Knepley   Notes:
825dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
826dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
8276858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
8286858538eSMatthew G. Knepley 
8290b3275a6SBarry Smith   The vector `c` should be destroyed by the caller.
8306858538eSMatthew G. Knepley 
831dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
8326858538eSMatthew G. Knepley @*/
833d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
834d71ae5a4SJacob Faibussowitsch {
8356858538eSMatthew G. Knepley   PetscFunctionBegin;
8366858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8376858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
8386858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
8396858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
8406858538eSMatthew G. Knepley   dm->coordinates[1].xl = c;
8416858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
8423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8436858538eSMatthew G. Knepley }
8446858538eSMatthew G. Knepley 
845d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
846d71ae5a4SJacob Faibussowitsch {
8476858538eSMatthew G. Knepley   PetscFunctionBegin;
8486858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8494f572ea9SToby Isaac   PetscAssertPointer(field, 2);
850213acdd3SPierre Jolivet   if (!dm->coordinates[0].field) PetscTryTypeMethod(dm, createcoordinatefield, &dm->coordinates[0].field);
8516858538eSMatthew G. Knepley   *field = dm->coordinates[0].field;
8523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8536858538eSMatthew G. Knepley }
8546858538eSMatthew G. Knepley 
855d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
856d71ae5a4SJacob Faibussowitsch {
8576858538eSMatthew G. Knepley   PetscFunctionBegin;
8586858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8596858538eSMatthew G. Knepley   if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2);
8606858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)field));
8616858538eSMatthew G. Knepley   PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
8626858538eSMatthew G. Knepley   dm->coordinates[0].field = field;
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8646858538eSMatthew G. Knepley }
8656858538eSMatthew G. Knepley 
8666c6a6b79SMatthew G. Knepley PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM dm, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[])
867d71ae5a4SJacob Faibussowitsch {
8686858538eSMatthew G. Knepley   Vec         coords = NULL;
8696858538eSMatthew G. Knepley   PetscReal   min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
8706858538eSMatthew G. Knepley   PetscReal   max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
8716858538eSMatthew G. Knepley   PetscInt    cdim, i, j;
8726c6a6b79SMatthew G. Knepley   PetscMPIInt size;
8736858538eSMatthew G. Knepley 
8746858538eSMatthew G. Knepley   PetscFunctionBegin;
8756c6a6b79SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
8766858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
8776c6a6b79SMatthew G. Knepley   if (size == 1) {
8786c6a6b79SMatthew G. Knepley     const PetscReal *L, *Lstart;
8796c6a6b79SMatthew G. Knepley 
8806c6a6b79SMatthew G. Knepley     PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
8816c6a6b79SMatthew G. Knepley     if (L) {
8826c6a6b79SMatthew G. Knepley       for (PetscInt d = 0; d < cdim; ++d)
8836c6a6b79SMatthew G. Knepley         if (L[d] > 0.0) {
8846c6a6b79SMatthew G. Knepley           min[d] = Lstart[d];
8856c6a6b79SMatthew G. Knepley           max[d] = Lstart[d] + L[d];
8866c6a6b79SMatthew G. Knepley         }
8876c6a6b79SMatthew G. Knepley     }
8886c6a6b79SMatthew G. Knepley   }
889bdf15c88SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
890bdf15c88SMatthew G. Knepley   if (coords) {
891bdf15c88SMatthew G. Knepley     const PetscScalar *local_coords;
892bdf15c88SMatthew G. Knepley     PetscInt           N, Ni;
893bdf15c88SMatthew G. Knepley 
894bdf15c88SMatthew G. Knepley     for (j = cdim; j < 3; ++j) {
895bdf15c88SMatthew G. Knepley       min[j] = 0;
896bdf15c88SMatthew G. Knepley       max[j] = 0;
897bdf15c88SMatthew G. Knepley     }
898bdf15c88SMatthew G. Knepley     PetscCall(VecGetArrayRead(coords, &local_coords));
899bdf15c88SMatthew G. Knepley     PetscCall(VecGetLocalSize(coords, &N));
900bdf15c88SMatthew G. Knepley     Ni = N / cdim;
901bdf15c88SMatthew G. Knepley     for (i = 0; i < Ni; ++i) {
902bdf15c88SMatthew G. Knepley       for (j = 0; j < cdim; ++j) {
903bdf15c88SMatthew G. Knepley         min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
904bdf15c88SMatthew G. Knepley         max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
905bdf15c88SMatthew G. Knepley       }
906bdf15c88SMatthew G. Knepley     }
907bdf15c88SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(coords, &local_coords));
908bdf15c88SMatthew G. Knepley     PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
9096858538eSMatthew G. Knepley     if (coords) {
9106858538eSMatthew G. Knepley       PetscCall(VecGetArrayRead(coords, &local_coords));
9116858538eSMatthew G. Knepley       PetscCall(VecGetLocalSize(coords, &N));
9126858538eSMatthew G. Knepley       Ni = N / cdim;
9136858538eSMatthew G. Knepley       for (i = 0; i < Ni; ++i) {
914bdf15c88SMatthew G. Knepley         for (j = 0; j < cdim; ++j) {
915bdf15c88SMatthew G. Knepley           min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
916bdf15c88SMatthew G. Knepley           max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
9176858538eSMatthew G. Knepley         }
9186858538eSMatthew G. Knepley       }
9196858538eSMatthew G. Knepley       PetscCall(VecRestoreArrayRead(coords, &local_coords));
920bdf15c88SMatthew G. Knepley     }
9216858538eSMatthew G. Knepley     if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
9226858538eSMatthew G. Knepley     if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
9236c6a6b79SMatthew G. Knepley   }
9246c6a6b79SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
9256c6a6b79SMatthew G. Knepley }
9266c6a6b79SMatthew G. Knepley 
9276c6a6b79SMatthew G. Knepley /*@
9286c6a6b79SMatthew G. Knepley   DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.
9296c6a6b79SMatthew G. Knepley 
9306c6a6b79SMatthew G. Knepley   Not Collective
9316c6a6b79SMatthew G. Knepley 
9326c6a6b79SMatthew G. Knepley   Input Parameter:
9336c6a6b79SMatthew G. Knepley . dm - the `DM`
9346c6a6b79SMatthew G. Knepley 
9356c6a6b79SMatthew G. Knepley   Output Parameters:
9366c6a6b79SMatthew G. Knepley + lmin - local minimum coordinates (length coord dim, optional)
9376c6a6b79SMatthew G. Knepley - lmax - local maximum coordinates (length coord dim, optional)
9386c6a6b79SMatthew G. Knepley 
9396c6a6b79SMatthew G. Knepley   Level: beginner
9406c6a6b79SMatthew G. Knepley 
9416c6a6b79SMatthew G. Knepley   Note:
9426c6a6b79SMatthew G. Knepley   If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.
9436c6a6b79SMatthew G. Knepley 
9446c6a6b79SMatthew G. Knepley .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
9456c6a6b79SMatthew G. Knepley @*/
9466c6a6b79SMatthew G. Knepley PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
9476c6a6b79SMatthew G. Knepley {
9486c6a6b79SMatthew G. Knepley   PetscFunctionBegin;
9496c6a6b79SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9506c6a6b79SMatthew G. Knepley   PetscUseTypeMethod(dm, getlocalboundingbox, lmin, lmax, NULL, NULL);
9513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9526858538eSMatthew G. Knepley }
9536858538eSMatthew G. Knepley 
9546858538eSMatthew G. Knepley /*@
955dce8aebaSBarry Smith   DMGetBoundingBox - Returns the global bounding box for the `DM`.
9566858538eSMatthew G. Knepley 
9576858538eSMatthew G. Knepley   Collective
9586858538eSMatthew G. Knepley 
9596858538eSMatthew G. Knepley   Input Parameter:
960dce8aebaSBarry Smith . dm - the `DM`
9616858538eSMatthew G. Knepley 
9626858538eSMatthew G. Knepley   Output Parameters:
9636858538eSMatthew G. Knepley + gmin - global minimum coordinates (length coord dim, optional)
96435cb6cd3SPierre Jolivet - gmax - global maximum coordinates (length coord dim, optional)
9656858538eSMatthew G. Knepley 
9666858538eSMatthew G. Knepley   Level: beginner
9676858538eSMatthew G. Knepley 
968dce8aebaSBarry Smith .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
9696858538eSMatthew G. Knepley @*/
970d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
971d71ae5a4SJacob Faibussowitsch {
9726858538eSMatthew G. Knepley   PetscReal        lmin[3], lmax[3];
9736c6a6b79SMatthew G. Knepley   const PetscReal *L, *Lstart;
9746858538eSMatthew G. Knepley   PetscInt         cdim;
9756858538eSMatthew G. Knepley 
9766858538eSMatthew G. Knepley   PetscFunctionBegin;
9776858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9786858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
9796858538eSMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
980e91c04dfSPierre Jolivet   if (gmin) PetscCallMPI(MPIU_Allreduce(lmin, gmin, cdim, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
981e91c04dfSPierre Jolivet   if (gmax) PetscCallMPI(MPIU_Allreduce(lmax, gmax, cdim, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
9826c6a6b79SMatthew G. Knepley   PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
9836c6a6b79SMatthew G. Knepley   if (L) {
9846c6a6b79SMatthew G. Knepley     for (PetscInt d = 0; d < cdim; ++d)
9856c6a6b79SMatthew G. Knepley       if (L[d] > 0.0) {
9866c6a6b79SMatthew G. Knepley         gmin[d] = Lstart[d];
9876c6a6b79SMatthew G. Knepley         gmax[d] = Lstart[d] + L[d];
9886c6a6b79SMatthew G. Knepley       }
9896c6a6b79SMatthew G. Knepley   }
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9916858538eSMatthew G. Knepley }
9926858538eSMatthew G. Knepley 
993e44f6aebSMatthew G. Knepley static PetscErrorCode DMCreateAffineCoordinates_Internal(DM dm)
99490612307SJed Brown {
995e44f6aebSMatthew G. Knepley   DM             cdm;
996e44f6aebSMatthew G. Knepley   PetscFE        feLinear;
997e44f6aebSMatthew G. Knepley   DMPolytopeType ct;
9981df12153SMatthew G. Knepley   PetscInt       dim, dE, height, cStart, cEnd, gct;
999e44f6aebSMatthew G. Knepley 
1000e44f6aebSMatthew G. Knepley   PetscFunctionBegin;
1001e44f6aebSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
1002e44f6aebSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
1003e44f6aebSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &dE));
10041df12153SMatthew G. Knepley   PetscCall(DMPlexGetVTKCellHeight(dm, &height));
10051df12153SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, height, &cStart, &cEnd));
1006e44f6aebSMatthew G. Knepley   if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1007e44f6aebSMatthew G. Knepley   else ct = DM_POLYTOPE_UNKNOWN;
1008e44f6aebSMatthew G. Knepley   gct = (PetscInt)ct;
1009462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &gct, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1010e44f6aebSMatthew G. Knepley   ct = (DMPolytopeType)gct;
1011e44f6aebSMatthew G. Knepley   // Work around current bug in PetscDualSpaceSetUp_Lagrange()
1012e44f6aebSMatthew G. Knepley   //   Can be seen in plex_tutorials-ex10_1
1013e44f6aebSMatthew G. Knepley   if (ct != DM_POLYTOPE_SEG_PRISM_TENSOR && ct != DM_POLYTOPE_TRI_PRISM_TENSOR && ct != DM_POLYTOPE_QUAD_PRISM_TENSOR) {
1014e44f6aebSMatthew G. Knepley     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, 1, -1, &feLinear));
1015e44f6aebSMatthew G. Knepley     PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)feLinear));
1016e44f6aebSMatthew G. Knepley     PetscCall(PetscFEDestroy(&feLinear));
1017e44f6aebSMatthew G. Knepley     PetscCall(DMCreateDS(cdm));
1018e44f6aebSMatthew G. Knepley   }
1019e44f6aebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1020e44f6aebSMatthew G. Knepley }
1021e44f6aebSMatthew G. Knepley 
1022e44f6aebSMatthew G. Knepley PetscErrorCode DMGetCoordinateDegree_Internal(DM dm, PetscInt *degree)
1023e44f6aebSMatthew G. Knepley {
1024e44f6aebSMatthew G. Knepley   DM           cdm;
1025e44f6aebSMatthew G. Knepley   PetscFE      fe;
1026e44f6aebSMatthew G. Knepley   PetscSpace   sp;
1027e44f6aebSMatthew G. Knepley   PetscClassId id;
1028e44f6aebSMatthew G. Knepley 
1029e44f6aebSMatthew G. Knepley   PetscFunctionBegin;
1030e44f6aebSMatthew G. Knepley   *degree = 1;
1031e44f6aebSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
1032e44f6aebSMatthew G. Knepley   PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
1033e44f6aebSMatthew G. Knepley   PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
1034e44f6aebSMatthew G. Knepley   if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
1035e44f6aebSMatthew G. Knepley   PetscCall(PetscFEGetBasisSpace(fe, &sp));
1036e44f6aebSMatthew G. Knepley   PetscCall(PetscSpaceGetDegree(sp, degree, NULL));
1037e44f6aebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
103890612307SJed Brown }
103990612307SJed Brown 
104090df3356SJames Wright static void evaluate_coordinates(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xnew[])
104190df3356SJames Wright {
104290df3356SJames Wright   for (PetscInt i = 0; i < dim; i++) xnew[i] = x[i];
104390df3356SJames Wright }
104490df3356SJames Wright 
10456858538eSMatthew G. Knepley /*@
1046e44f6aebSMatthew G. Knepley   DMSetCoordinateDisc - Set a coordinate space
10476858538eSMatthew G. Knepley 
10486858538eSMatthew G. Knepley   Input Parameters:
1049dce8aebaSBarry Smith + dm      - The `DM` object
10504688473bSSatish Balay . disc    - The new coordinate discretization or `NULL` to ensure a coordinate discretization exists
1051e44f6aebSMatthew G. Knepley - project - Project coordinates to new discretization
10526858538eSMatthew G. Knepley 
10536858538eSMatthew G. Knepley   Level: intermediate
10546858538eSMatthew G. Knepley 
1055dce8aebaSBarry Smith   Notes:
1056e44f6aebSMatthew G. Knepley   A `PetscFE` defines an approximation space using a `PetscSpace`, which represents the basis functions, and a `PetscDualSpace`, which defines the interpolation operation in the space.
1057dce8aebaSBarry Smith 
1058dce8aebaSBarry Smith   This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
1059e44f6aebSMatthew G. Knepley   The coordinate projection is done on the continuous coordinates, but the discontinuous coordinates are not updated.
1060dce8aebaSBarry Smith 
10610b3275a6SBarry Smith   Developer Note:
1062dce8aebaSBarry Smith   With more effort, we could directly project the discontinuous coordinates also.
1063dce8aebaSBarry Smith 
1064dce8aebaSBarry Smith .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
10656858538eSMatthew G. Knepley @*/
1066e44f6aebSMatthew G. Knepley PetscErrorCode DMSetCoordinateDisc(DM dm, PetscFE disc, PetscBool project)
1067d71ae5a4SJacob Faibussowitsch {
1068e44f6aebSMatthew G. Knepley   DM           cdmOld, cdmNew;
10696858538eSMatthew G. Knepley   PetscFE      discOld;
10706858538eSMatthew G. Knepley   PetscClassId classid;
10716858538eSMatthew G. Knepley   PetscBool    same_space = PETSC_TRUE;
1072dd4c3f67SMatthew G. Knepley   const char  *prefix;
10736858538eSMatthew G. Knepley 
10746858538eSMatthew G. Knepley   PetscFunctionBegin;
10756858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10766858538eSMatthew G. Knepley   if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2);
10776858538eSMatthew G. Knepley 
10786858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdmOld));
10796858538eSMatthew G. Knepley   /* Check current discretization is compatible */
10806858538eSMatthew G. Knepley   PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
10816858538eSMatthew G. Knepley   PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
10826858538eSMatthew G. Knepley   if (classid != PETSCFE_CLASSID) {
10836858538eSMatthew G. Knepley     if (classid == PETSC_CONTAINER_CLASSID) {
1084e44f6aebSMatthew G. Knepley       PetscCall(DMCreateAffineCoordinates_Internal(dm));
10856858538eSMatthew G. Knepley       PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
10866858538eSMatthew G. Knepley     } else {
10876858538eSMatthew G. Knepley       const char *discname;
10886858538eSMatthew G. Knepley 
10896858538eSMatthew G. Knepley       PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
10906858538eSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
10916858538eSMatthew G. Knepley     }
10926858538eSMatthew G. Knepley   }
1093e44f6aebSMatthew G. Knepley   // Linear space has been created by now
10943ba16761SJacob Faibussowitsch   if (!disc) PetscFunctionReturn(PETSC_SUCCESS);
1095e44f6aebSMatthew G. Knepley   // Check if the new space is the same as the old modulo quadrature
1096e44f6aebSMatthew G. Knepley   {
10976858538eSMatthew G. Knepley     PetscDualSpace dsOld, ds;
10986858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
10996858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(disc, &ds));
11006858538eSMatthew G. Knepley     PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
11016858538eSMatthew G. Knepley   }
1102e44f6aebSMatthew G. Knepley   // Make a fresh clone of the coordinate DM
11036858538eSMatthew G. Knepley   PetscCall(DMClone(cdmOld, &cdmNew));
1104dd4c3f67SMatthew G. Knepley   cdmNew->cloneOpts = PETSC_TRUE;
1105dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix));
1106dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix));
11076858538eSMatthew G. Knepley   PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
11086858538eSMatthew G. Knepley   PetscCall(DMCreateDS(cdmNew));
1109e44f6aebSMatthew G. Knepley   {
1110e44f6aebSMatthew G. Knepley     PetscDS ds, nds;
1111e44f6aebSMatthew G. Knepley 
1112e44f6aebSMatthew G. Knepley     PetscCall(DMGetDS(cdmOld, &ds));
1113e44f6aebSMatthew G. Knepley     PetscCall(DMGetDS(cdmNew, &nds));
1114e44f6aebSMatthew G. Knepley     PetscCall(PetscDSCopyConstants(ds, nds));
1115e44f6aebSMatthew G. Knepley   }
111681316d9bSJames Wright   if (cdmOld->periodic.setup) {
111781316d9bSJames Wright     PetscSF dummy;
111881316d9bSJames Wright     // Force IsoperiodicPointSF to be built, required for periodic coordinate setup
111981316d9bSJames Wright     PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &dummy));
112081316d9bSJames Wright     cdmNew->periodic.setup = cdmOld->periodic.setup;
112181316d9bSJames Wright     PetscCall(cdmNew->periodic.setup(cdmNew));
112281316d9bSJames Wright   }
1123dd4c3f67SMatthew G. Knepley   if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew));
1124e44f6aebSMatthew G. Knepley   if (project) {
112553b735f8SJames Wright     Vec      coordsOld, coordsNew;
112681316d9bSJames Wright     PetscInt num_face_sfs = 0;
112753b735f8SJames Wright 
112881316d9bSJames Wright     PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, NULL));
112981316d9bSJames Wright     if (num_face_sfs) { // Isoperiodicity requires projecting the local coordinates
113053b735f8SJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &coordsOld));
113153b735f8SJames Wright       PetscCall(DMCreateLocalVector(cdmNew, &coordsNew));
113253b735f8SJames Wright       PetscCall(PetscObjectSetName((PetscObject)coordsNew, "coordinates"));
1133f0ac6e35SMatthew G. Knepley       if (same_space) {
1134f0ac6e35SMatthew G. Knepley         // Need to copy so that the new vector has the right dm
1135f0ac6e35SMatthew G. Knepley         PetscCall(VecCopy(coordsOld, coordsNew));
1136e44f6aebSMatthew G. Knepley       } else {
113790df3356SJames Wright         void (*funcs[])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {evaluate_coordinates};
113890df3356SJames Wright 
113990df3356SJames Wright         // We can't call DMProjectField directly because it depends on KSP for DMGlobalToLocalSolve(), but we can use the core strategy
114090df3356SJames Wright         PetscCall(DMSetCoordinateDM(cdmNew, cdmOld));
114190df3356SJames Wright         // See DMPlexRemapGeometry() for a similar pattern handling the coordinate field
114290df3356SJames Wright         DMField cf;
114390df3356SJames Wright         PetscCall(DMGetCoordinateField(dm, &cf));
114490df3356SJames Wright         cdmNew->coordinates[0].field = cf;
114553b735f8SJames Wright         PetscCall(DMProjectFieldLocal(cdmNew, 0.0, NULL, funcs, INSERT_VALUES, coordsNew));
114690df3356SJames Wright         cdmNew->coordinates[0].field = NULL;
114790df3356SJames Wright         PetscCall(DMSetCoordinateDM(cdmNew, NULL));
114853b735f8SJames Wright       }
114953b735f8SJames Wright       PetscCall(DMSetCoordinatesLocal(dm, coordsNew));
115053b735f8SJames Wright       PetscCall(VecDestroy(&coordsNew));
115153b735f8SJames Wright     } else {
115253b735f8SJames Wright       PetscCall(DMGetCoordinates(dm, &coordsOld));
115353b735f8SJames Wright       PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
115453b735f8SJames Wright       if (same_space) {
115553b735f8SJames Wright         // Need to copy so that the new vector has the right dm
115653b735f8SJames Wright         PetscCall(VecCopy(coordsOld, coordsNew));
115790df3356SJames Wright       } else {
1158e44f6aebSMatthew G. Knepley         Mat In;
1159e44f6aebSMatthew G. Knepley 
1160e44f6aebSMatthew G. Knepley         PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &In, NULL));
1161e44f6aebSMatthew G. Knepley         PetscCall(MatMult(In, coordsOld, coordsNew));
1162e44f6aebSMatthew G. Knepley         PetscCall(MatDestroy(&In));
1163e44f6aebSMatthew G. Knepley       }
1164e44f6aebSMatthew G. Knepley       PetscCall(DMSetCoordinates(dm, coordsNew));
1165e44f6aebSMatthew G. Knepley       PetscCall(VecDestroy(&coordsNew));
11666858538eSMatthew G. Knepley     }
116753b735f8SJames Wright   }
11686858538eSMatthew G. Knepley   /* Set new coordinate structures */
11696858538eSMatthew G. Knepley   PetscCall(DMSetCoordinateField(dm, NULL));
11706858538eSMatthew G. Knepley   PetscCall(DMSetCoordinateDM(dm, cdmNew));
11716858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cdmNew));
11723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11736858538eSMatthew G. Knepley }
11746858538eSMatthew G. Knepley 
11756858538eSMatthew G. Knepley /*@
11760b3275a6SBarry Smith   DMLocatePoints - Locate the points in `v` in the mesh and return a `PetscSF` of the containing cells
11776858538eSMatthew G. Knepley 
117820f4b53cSBarry Smith   Collective
11796858538eSMatthew G. Knepley 
11806858538eSMatthew G. Knepley   Input Parameters:
1181dce8aebaSBarry Smith + dm    - The `DM`
1182dce8aebaSBarry Smith - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`
11836858538eSMatthew G. Knepley 
11846858538eSMatthew G. Knepley   Input/Output Parameters:
1185dce8aebaSBarry Smith + v      - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
118620f4b53cSBarry Smith - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point;
11870b3275a6SBarry Smith            on output, the `PetscSF` containing the MPI ranks and local indices of the containing points
11886858538eSMatthew G. Knepley 
11896858538eSMatthew G. Knepley   Level: developer
11906858538eSMatthew G. Knepley 
11916858538eSMatthew G. Knepley   Notes:
11920b3275a6SBarry Smith   To do a search of the local cells of the mesh, `v` should have `PETSC_COMM_SELF` as its communicator.
11930b3275a6SBarry Smith   To do a search of all the cells in the distributed mesh, `v` should have the same MPI communicator as `dm`.
11946858538eSMatthew G. Knepley 
1195d8206211SMatthew G. Knepley   Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.
1196d8206211SMatthew G. Knepley 
119720f4b53cSBarry Smith   If *cellSF is `NULL` on input, a `PetscSF` will be created.
119820f4b53cSBarry Smith   If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.
11996858538eSMatthew G. Knepley 
12006858538eSMatthew G. Knepley   An array that maps each point to its containing cell can be obtained with
1201dce8aebaSBarry Smith .vb
1202dce8aebaSBarry Smith     const PetscSFNode *cells;
1203dce8aebaSBarry Smith     PetscInt           nFound;
1204dce8aebaSBarry Smith     const PetscInt    *found;
12056858538eSMatthew G. Knepley 
1206dce8aebaSBarry Smith     PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1207dce8aebaSBarry Smith .ve
12086858538eSMatthew G. Knepley 
12090b3275a6SBarry Smith   Where cells[i].rank is the MPI rank of the process owning the cell containing point found[i] (or i if found == NULL), and cells[i].index is
12100b3275a6SBarry Smith   the index of the cell in its MPI process' local numbering. This rank is in the communicator for `v`, so if `v` is on `PETSC_COMM_SELF` then the rank will always be 0.
12116858538eSMatthew G. Knepley 
1212dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
12136858538eSMatthew G. Knepley @*/
1214d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1215d71ae5a4SJacob Faibussowitsch {
12166858538eSMatthew G. Knepley   PetscFunctionBegin;
12176858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12186858538eSMatthew G. Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 2);
12194f572ea9SToby Isaac   PetscAssertPointer(cellSF, 4);
12206858538eSMatthew G. Knepley   if (*cellSF) {
12216858538eSMatthew G. Knepley     PetscMPIInt result;
12226858538eSMatthew G. Knepley 
12236858538eSMatthew G. Knepley     PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4);
12246858538eSMatthew G. Knepley     PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
12256858538eSMatthew G. Knepley     PetscCheck(result == MPI_IDENT || result == MPI_CONGRUENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "cellSF must have a communicator congruent to v's");
12266858538eSMatthew G. Knepley   } else {
12276858538eSMatthew G. Knepley     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
12286858538eSMatthew G. Knepley   }
12296858538eSMatthew G. Knepley   PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1230dbbe0bcdSBarry Smith   PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
12316858538eSMatthew G. Knepley   PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
12323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12336858538eSMatthew G. Knepley }
1234