xref: /petsc/src/dm/interface/dmcoordinates.c (revision dce8aeba1c9b69b19f651c53d8a6b674bd7e9cbd)
16858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
26858538eSMatthew G. Knepley 
36858538eSMatthew G. Knepley #include <petscdmplex.h> /* For DMProjectCoordinates() */
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;
116858538eSMatthew G. Knepley   PetscFunctionBegin;
126858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
136858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dmc, &dmc_coord));
146858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
156858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dmc, &ccoords));
166858538eSMatthew G. Knepley   if (coords && !ccoords) {
176858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords));
186858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
196858538eSMatthew G. Knepley     PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject));
206858538eSMatthew G. Knepley     PetscCall(MatRestrict(inject, coords, ccoords));
216858538eSMatthew G. Knepley     PetscCall(MatDestroy(&inject));
226858538eSMatthew G. Knepley     PetscCall(DMSetCoordinates(dmc, ccoords));
236858538eSMatthew G. Knepley     PetscCall(VecDestroy(&ccoords));
246858538eSMatthew G. Knepley   }
256858538eSMatthew G. Knepley   PetscFunctionReturn(0);
266858538eSMatthew G. Knepley }
276858538eSMatthew G. Knepley 
28d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, void *ctx)
29d71ae5a4SJacob Faibussowitsch {
306858538eSMatthew G. Knepley   DM          dm_coord, subdm_coord;
316858538eSMatthew G. Knepley   Vec         coords, ccoords, clcoords;
326858538eSMatthew G. Knepley   VecScatter *scat_i, *scat_g;
336858538eSMatthew G. Knepley   PetscFunctionBegin;
346858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
356858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(subdm, &subdm_coord));
366858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
376858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(subdm, &ccoords));
386858538eSMatthew G. Knepley   if (coords && !ccoords) {
396858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords));
406858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
416858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(subdm_coord, &clcoords));
426858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates"));
436858538eSMatthew G. Knepley     PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g));
446858538eSMatthew G. Knepley     PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
456858538eSMatthew G. Knepley     PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
466858538eSMatthew G. Knepley     PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
476858538eSMatthew G. Knepley     PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
486858538eSMatthew G. Knepley     PetscCall(DMSetCoordinates(subdm, ccoords));
496858538eSMatthew G. Knepley     PetscCall(DMSetCoordinatesLocal(subdm, clcoords));
506858538eSMatthew G. Knepley     PetscCall(VecScatterDestroy(&scat_i[0]));
516858538eSMatthew G. Knepley     PetscCall(VecScatterDestroy(&scat_g[0]));
526858538eSMatthew G. Knepley     PetscCall(VecDestroy(&ccoords));
536858538eSMatthew G. Knepley     PetscCall(VecDestroy(&clcoords));
546858538eSMatthew G. Knepley     PetscCall(PetscFree(scat_i));
556858538eSMatthew G. Knepley     PetscCall(PetscFree(scat_g));
566858538eSMatthew G. Knepley   }
576858538eSMatthew G. Knepley   PetscFunctionReturn(0);
586858538eSMatthew G. Knepley }
596858538eSMatthew G. Knepley 
606858538eSMatthew G. Knepley /*@
61*dce8aebaSBarry Smith   DMGetCoordinateDM - Gets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
626858538eSMatthew G. Knepley 
636858538eSMatthew G. Knepley   Collective on dm
646858538eSMatthew G. Knepley 
656858538eSMatthew G. Knepley   Input Parameter:
66*dce8aebaSBarry Smith . dm - the `DM`
676858538eSMatthew G. Knepley 
686858538eSMatthew G. Knepley   Output Parameter:
69*dce8aebaSBarry Smith . cdm - coordinate `DM`
706858538eSMatthew G. Knepley 
716858538eSMatthew G. Knepley   Level: intermediate
726858538eSMatthew G. Knepley 
73*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGSetCellCoordinateDM()`,
74*dce8aebaSBarry Smith           `DMGSetCellCoordinateDM()`
756858538eSMatthew G. Knepley @*/
76d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
77d71ae5a4SJacob Faibussowitsch {
786858538eSMatthew G. Knepley   PetscFunctionBegin;
796858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
806858538eSMatthew G. Knepley   PetscValidPointer(cdm, 2);
816858538eSMatthew G. Knepley   if (!dm->coordinates[0].dm) {
826858538eSMatthew G. Knepley     DM cdm;
836858538eSMatthew G. Knepley 
84dbbe0bcdSBarry Smith     PetscUseTypeMethod(dm, createcoordinatedm, &cdm);
856858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM"));
866858538eSMatthew G. Knepley     /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup
876858538eSMatthew G. Knepley      * until the call to CreateCoordinateDM) */
886858538eSMatthew G. Knepley     PetscCall(DMDestroy(&dm->coordinates[0].dm));
896858538eSMatthew G. Knepley     dm->coordinates[0].dm = cdm;
906858538eSMatthew G. Knepley   }
916858538eSMatthew G. Knepley   *cdm = dm->coordinates[0].dm;
926858538eSMatthew G. Knepley   PetscFunctionReturn(0);
936858538eSMatthew G. Knepley }
946858538eSMatthew G. Knepley 
956858538eSMatthew G. Knepley /*@
96*dce8aebaSBarry Smith   DMSetCoordinateDM - Sets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
976858538eSMatthew G. Knepley 
986858538eSMatthew G. Knepley   Logically Collective on dm
996858538eSMatthew G. Knepley 
1006858538eSMatthew G. Knepley   Input Parameters:
101*dce8aebaSBarry Smith + dm - the `DM`
102*dce8aebaSBarry Smith - cdm - coordinate `DM`
1036858538eSMatthew G. Knepley 
1046858538eSMatthew G. Knepley   Level: intermediate
1056858538eSMatthew G. Knepley 
106*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMGetCellCoordinateDM()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`,
107*dce8aebaSBarry Smith           `DMGSetCellCoordinateDM()`
1086858538eSMatthew G. Knepley @*/
109d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm)
110d71ae5a4SJacob Faibussowitsch {
1116858538eSMatthew G. Knepley   PetscFunctionBegin;
1126858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1136858538eSMatthew G. Knepley   PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
1146858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)cdm));
1156858538eSMatthew G. Knepley   PetscCall(DMDestroy(&dm->coordinates[0].dm));
1166858538eSMatthew G. Knepley   dm->coordinates[0].dm = cdm;
1176858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1186858538eSMatthew G. Knepley }
1196858538eSMatthew G. Knepley 
1206858538eSMatthew G. Knepley /*@
121*dce8aebaSBarry Smith   DMGetCellCoordinateDM - Gets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1226858538eSMatthew G. Knepley 
1236858538eSMatthew G. Knepley   Collective on dm
1246858538eSMatthew G. Knepley 
1256858538eSMatthew G. Knepley   Input Parameter:
126*dce8aebaSBarry Smith . dm - the `DM`
1276858538eSMatthew G. Knepley 
1286858538eSMatthew G. Knepley   Output Parameter:
129*dce8aebaSBarry Smith . cdm - cellwise coordinate `DM`, or NULL if they are not defined
1306858538eSMatthew G. Knepley 
1316858538eSMatthew G. Knepley   Level: intermediate
1326858538eSMatthew G. Knepley 
133*dce8aebaSBarry Smith   Note:
134*dce8aebaSBarry Smith   Call `DMLocalizeCoordinates()` to automatically create cellwise coordinates for periodic geometries.
135*dce8aebaSBarry Smith 
136*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
137*dce8aebaSBarry Smith           `DMLocalizeCoordinates()`, `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
1386858538eSMatthew G. Knepley @*/
139d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm)
140d71ae5a4SJacob Faibussowitsch {
1416858538eSMatthew G. Knepley   PetscFunctionBegin;
1426858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1436858538eSMatthew G. Knepley   PetscValidPointer(cdm, 2);
1446858538eSMatthew G. Knepley   *cdm = dm->coordinates[1].dm;
1456858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1466858538eSMatthew G. Knepley }
1476858538eSMatthew G. Knepley 
1486858538eSMatthew G. Knepley /*@
149*dce8aebaSBarry Smith   DMSetCellCoordinateDM - Sets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1506858538eSMatthew G. Knepley 
1516858538eSMatthew G. Knepley   Logically Collective on dm
1526858538eSMatthew G. Knepley 
1536858538eSMatthew G. Knepley   Input Parameters:
154*dce8aebaSBarry Smith + dm - the `DM`
155*dce8aebaSBarry Smith - cdm - cellwise coordinate `DM`
1566858538eSMatthew G. Knepley 
1576858538eSMatthew G. Knepley   Level: intermediate
1586858538eSMatthew G. Knepley 
159*dce8aebaSBarry Smith   Note:
160*dce8aebaSBarry Smith   As opposed to `DMSetCoordinateDM()` these coordinates are useful for discontinous Galerkin methods since they support coordinate fields that are discontinuous at cell boundaries.
161*dce8aebaSBarry Smith 
162*dce8aebaSBarry Smith .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
163*dce8aebaSBarry Smith           `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
1646858538eSMatthew G. Knepley @*/
165d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm)
166d71ae5a4SJacob Faibussowitsch {
1676858538eSMatthew G. Knepley   PetscInt dim;
1686858538eSMatthew G. Knepley 
1696858538eSMatthew G. Knepley   PetscFunctionBegin;
1706858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1716858538eSMatthew G. Knepley   if (cdm) {
1726858538eSMatthew G. Knepley     PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
1736858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &dim));
1746858538eSMatthew G. Knepley     dm->coordinates[1].dim = dim;
1756858538eSMatthew G. Knepley   }
1766858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)cdm));
1776858538eSMatthew G. Knepley   PetscCall(DMDestroy(&dm->coordinates[1].dm));
1786858538eSMatthew G. Knepley   dm->coordinates[1].dm = cdm;
1796858538eSMatthew G. Knepley   PetscFunctionReturn(0);
1806858538eSMatthew G. Knepley }
1816858538eSMatthew G. Knepley 
1826858538eSMatthew G. Knepley /*@
183*dce8aebaSBarry Smith   DMGetCoordinateDim - Retrieve the dimension of embedding space for coordinate values. For example a mesh on the surface of a sphere would have a 3 dimensional embedding space
1846858538eSMatthew G. Knepley 
1856858538eSMatthew G. Knepley   Not Collective
1866858538eSMatthew G. Knepley 
1876858538eSMatthew G. Knepley   Input Parameter:
188*dce8aebaSBarry Smith . dm - The `DM` object
1896858538eSMatthew G. Knepley 
1906858538eSMatthew G. Knepley   Output Parameter:
1916858538eSMatthew G. Knepley . dim - The embedding dimension
1926858538eSMatthew G. Knepley 
1936858538eSMatthew G. Knepley   Level: intermediate
1946858538eSMatthew G. Knepley 
195*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
1966858538eSMatthew G. Knepley @*/
197d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim)
198d71ae5a4SJacob Faibussowitsch {
1996858538eSMatthew G. Knepley   PetscFunctionBegin;
2006858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2016858538eSMatthew G. Knepley   PetscValidIntPointer(dim, 2);
2026858538eSMatthew G. Knepley   if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim;
2036858538eSMatthew G. Knepley   *dim = dm->coordinates[0].dim;
2046858538eSMatthew G. Knepley   PetscFunctionReturn(0);
2056858538eSMatthew G. Knepley }
2066858538eSMatthew G. Knepley 
2076858538eSMatthew G. Knepley /*@
2086858538eSMatthew G. Knepley   DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values.
2096858538eSMatthew G. Knepley 
2106858538eSMatthew G. Knepley   Not Collective
2116858538eSMatthew G. Knepley 
2126858538eSMatthew G. Knepley   Input Parameters:
213*dce8aebaSBarry Smith + dm  - The `DM` object
2146858538eSMatthew G. Knepley - dim - The embedding dimension
2156858538eSMatthew G. Knepley 
2166858538eSMatthew G. Knepley   Level: intermediate
2176858538eSMatthew G. Knepley 
218*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2196858538eSMatthew G. Knepley @*/
220d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim)
221d71ae5a4SJacob Faibussowitsch {
2226858538eSMatthew G. Knepley   PetscDS  ds;
2236858538eSMatthew G. Knepley   PetscInt Nds, n;
2246858538eSMatthew G. Knepley 
2256858538eSMatthew G. Knepley   PetscFunctionBegin;
2266858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2276858538eSMatthew G. Knepley   dm->coordinates[0].dim = dim;
2286858538eSMatthew G. Knepley   if (dm->dim >= 0) {
2296858538eSMatthew G. Knepley     PetscCall(DMGetNumDS(dm, &Nds));
2306858538eSMatthew G. Knepley     for (n = 0; n < Nds; ++n) {
2316858538eSMatthew G. Knepley       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds));
2326858538eSMatthew G. Knepley       PetscCall(PetscDSSetCoordinateDimension(ds, dim));
2336858538eSMatthew G. Knepley     }
2346858538eSMatthew G. Knepley   }
2356858538eSMatthew G. Knepley   PetscFunctionReturn(0);
2366858538eSMatthew G. Knepley }
2376858538eSMatthew G. Knepley 
2386858538eSMatthew G. Knepley /*@
2396858538eSMatthew G. Knepley   DMGetCoordinateSection - Retrieve the layout of coordinate values over the mesh.
2406858538eSMatthew G. Knepley 
2416858538eSMatthew G. Knepley   Collective on dm
2426858538eSMatthew G. Knepley 
2436858538eSMatthew G. Knepley   Input Parameter:
244*dce8aebaSBarry Smith . dm - The `DM` object
2456858538eSMatthew G. Knepley 
2466858538eSMatthew G. Knepley   Output Parameter:
247*dce8aebaSBarry Smith . section - The `PetscSection` object
2486858538eSMatthew G. Knepley 
2496858538eSMatthew G. Knepley   Level: intermediate
2506858538eSMatthew G. Knepley 
251*dce8aebaSBarry Smith   Note:
252*dce8aebaSBarry Smith   This just retrieves the local section from the coordinate `DM`. In other words,
253*dce8aebaSBarry Smith .vb
254*dce8aebaSBarry Smith   DMGetCoordinateDM(dm, &cdm);
255*dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
256*dce8aebaSBarry Smith .ve
257*dce8aebaSBarry Smith 
2586858538eSMatthew G. Knepley .seealso: `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2596858538eSMatthew G. Knepley @*/
260d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section)
261d71ae5a4SJacob Faibussowitsch {
2626858538eSMatthew G. Knepley   DM cdm;
2636858538eSMatthew G. Knepley 
2646858538eSMatthew G. Knepley   PetscFunctionBegin;
2656858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2666858538eSMatthew G. Knepley   PetscValidPointer(section, 2);
2676858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2686858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, section));
2696858538eSMatthew G. Knepley   PetscFunctionReturn(0);
2706858538eSMatthew G. Knepley }
2716858538eSMatthew G. Knepley 
2726858538eSMatthew G. Knepley /*@
2736858538eSMatthew G. Knepley   DMSetCoordinateSection - Set the layout of coordinate values over the mesh.
2746858538eSMatthew G. Knepley 
2756858538eSMatthew G. Knepley   Not Collective
2766858538eSMatthew G. Knepley 
2776858538eSMatthew G. Knepley   Input Parameters:
278*dce8aebaSBarry Smith + dm      - The `DM` object
279*dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
280*dce8aebaSBarry Smith - section - The `PetscSection` object
2816858538eSMatthew G. Knepley 
2826858538eSMatthew G. Knepley   Level: intermediate
2836858538eSMatthew G. Knepley 
284*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2856858538eSMatthew G. Knepley @*/
286d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section)
287d71ae5a4SJacob Faibussowitsch {
2886858538eSMatthew G. Knepley   DM cdm;
2896858538eSMatthew G. Knepley 
2906858538eSMatthew G. Knepley   PetscFunctionBegin;
2916858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2926858538eSMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
2936858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2946858538eSMatthew G. Knepley   PetscCall(DMSetLocalSection(cdm, section));
2956858538eSMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
2966858538eSMatthew G. Knepley     PetscInt d = PETSC_DEFAULT;
2976858538eSMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
2986858538eSMatthew G. Knepley 
2996858538eSMatthew G. Knepley     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3006858538eSMatthew G. Knepley     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
3016858538eSMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
3026858538eSMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
3036858538eSMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
3046858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(section, v, &dd));
3059371c9d4SSatish Balay       if (dd) {
3069371c9d4SSatish Balay         d = dd;
3079371c9d4SSatish Balay         break;
3089371c9d4SSatish Balay       }
3096858538eSMatthew G. Knepley     }
3106858538eSMatthew G. Knepley     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
3116858538eSMatthew G. Knepley   }
3126858538eSMatthew G. Knepley   PetscFunctionReturn(0);
3136858538eSMatthew G. Knepley }
3146858538eSMatthew G. Knepley 
3156858538eSMatthew G. Knepley /*@
3166858538eSMatthew G. Knepley   DMGetCellCoordinateSection - Retrieve the layout of cellwise coordinate values over the mesh.
3176858538eSMatthew G. Knepley 
3186858538eSMatthew G. Knepley   Collective on dm
3196858538eSMatthew G. Knepley 
3206858538eSMatthew G. Knepley   Input Parameter:
321*dce8aebaSBarry Smith . dm - The `DM` object
3226858538eSMatthew G. Knepley 
3236858538eSMatthew G. Knepley   Output Parameter:
324*dce8aebaSBarry Smith . section - The `PetscSection` object, or NULL if no cellwise coordinates are defined
3256858538eSMatthew G. Knepley 
3266858538eSMatthew G. Knepley   Level: intermediate
3276858538eSMatthew G. Knepley 
328*dce8aebaSBarry Smith   Note:
329*dce8aebaSBarry Smith   This just retrieves the local section from the cell coordinate `DM`. In other words,
330*dce8aebaSBarry Smith .vb
331*dce8aebaSBarry Smith   DMGetCellCoordinateDM(dm, &cdm);
332*dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
333*dce8aebaSBarry Smith .ve
334*dce8aebaSBarry Smith 
335*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
3366858538eSMatthew G. Knepley @*/
337d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section)
338d71ae5a4SJacob Faibussowitsch {
3396858538eSMatthew G. Knepley   DM cdm;
3406858538eSMatthew G. Knepley 
3416858538eSMatthew G. Knepley   PetscFunctionBegin;
3426858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3436858538eSMatthew G. Knepley   PetscValidPointer(section, 2);
3446858538eSMatthew G. Knepley   *section = NULL;
3456858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3466858538eSMatthew G. Knepley   if (cdm) PetscCall(DMGetLocalSection(cdm, section));
3476858538eSMatthew G. Knepley   PetscFunctionReturn(0);
3486858538eSMatthew G. Knepley }
3496858538eSMatthew G. Knepley 
3506858538eSMatthew G. Knepley /*@
3516858538eSMatthew G. Knepley   DMSetCellCoordinateSection - Set the layout of cellwise coordinate values over the mesh.
3526858538eSMatthew G. Knepley 
3536858538eSMatthew G. Knepley   Not Collective
3546858538eSMatthew G. Knepley 
3556858538eSMatthew G. Knepley   Input Parameters:
356*dce8aebaSBarry Smith + dm      - The `DM` object
357*dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
358*dce8aebaSBarry Smith - section - The `PetscSection` object for a cellwise layout
3596858538eSMatthew G. Knepley 
3606858538eSMatthew G. Knepley   Level: intermediate
3616858538eSMatthew G. Knepley 
362*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
3636858538eSMatthew G. Knepley @*/
364d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section)
365d71ae5a4SJacob Faibussowitsch {
3666858538eSMatthew G. Knepley   DM cdm;
3676858538eSMatthew G. Knepley 
3686858538eSMatthew G. Knepley   PetscFunctionBegin;
3696858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3706858538eSMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
3716858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3726858538eSMatthew G. Knepley   PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates");
3736858538eSMatthew G. Knepley   PetscCall(DMSetLocalSection(cdm, section));
3746858538eSMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
3756858538eSMatthew G. Knepley     PetscInt d = PETSC_DEFAULT;
3766858538eSMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
3776858538eSMatthew G. Knepley 
3786858538eSMatthew G. Knepley     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3796858538eSMatthew G. Knepley     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
3806858538eSMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
3816858538eSMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
3826858538eSMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
3836858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(section, v, &dd));
3849371c9d4SSatish Balay       if (dd) {
3859371c9d4SSatish Balay         d = dd;
3869371c9d4SSatish Balay         break;
3879371c9d4SSatish Balay       }
3886858538eSMatthew G. Knepley     }
3896858538eSMatthew G. Knepley     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
3906858538eSMatthew G. Knepley   }
3916858538eSMatthew G. Knepley   PetscFunctionReturn(0);
3926858538eSMatthew G. Knepley }
3936858538eSMatthew G. Knepley 
3946858538eSMatthew G. Knepley /*@
395*dce8aebaSBarry Smith   DMGetCoordinates - Gets a global vector with the coordinates associated with the `DM`.
3966858538eSMatthew G. Knepley 
3976858538eSMatthew G. Knepley   Collective on dm
3986858538eSMatthew G. Knepley 
3996858538eSMatthew G. Knepley   Input Parameter:
400*dce8aebaSBarry Smith . dm - the `DM`
4016858538eSMatthew G. Knepley 
4026858538eSMatthew G. Knepley   Output Parameter:
4036858538eSMatthew G. Knepley . c - global coordinate vector
4046858538eSMatthew G. Knepley 
405*dce8aebaSBarry Smith   Level: intermediate
406*dce8aebaSBarry Smith 
407*dce8aebaSBarry Smith   Notes:
408*dce8aebaSBarry Smith   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
4096858538eSMatthew G. Knepley   destroyed the array will no longer be valid.
4106858538eSMatthew G. Knepley 
4116858538eSMatthew G. Knepley   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
4126858538eSMatthew G. Knepley 
413*dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
4146858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
4156858538eSMatthew G. Knepley 
416*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
4176858538eSMatthew G. Knepley @*/
418d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
419d71ae5a4SJacob Faibussowitsch {
4206858538eSMatthew G. Knepley   PetscFunctionBegin;
4216858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4226858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
4236858538eSMatthew G. Knepley   if (!dm->coordinates[0].x && dm->coordinates[0].xl) {
4246858538eSMatthew G. Knepley     DM cdm = NULL;
4256858538eSMatthew G. Knepley 
4266858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
4276858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x));
4286858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates"));
4296858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
4306858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
4316858538eSMatthew G. Knepley   }
4326858538eSMatthew G. Knepley   *c = dm->coordinates[0].x;
4336858538eSMatthew G. Knepley   PetscFunctionReturn(0);
4346858538eSMatthew G. Knepley }
4356858538eSMatthew G. Knepley 
4366858538eSMatthew G. Knepley /*@
437*dce8aebaSBarry Smith   DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates
4386858538eSMatthew G. Knepley 
4396858538eSMatthew G. Knepley   Collective on dm
4406858538eSMatthew G. Knepley 
4416858538eSMatthew G. Knepley   Input Parameters:
442*dce8aebaSBarry Smith + dm - the `DM`
4436858538eSMatthew G. Knepley - c - coordinate vector
4446858538eSMatthew G. Knepley 
4456858538eSMatthew G. Knepley   Level: intermediate
4466858538eSMatthew G. Knepley 
447*dce8aebaSBarry Smith   Notes:
448*dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
449*dce8aebaSBarry Smith 
450*dce8aebaSBarry Smith   The vector c can be destroyed after the call
451*dce8aebaSBarry Smith 
452*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
4536858538eSMatthew G. Knepley @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinates(DM dm, Vec c)
455d71ae5a4SJacob Faibussowitsch {
4566858538eSMatthew G. Knepley   PetscFunctionBegin;
4576858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4586858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
4596858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
4606858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
4616858538eSMatthew G. Knepley   dm->coordinates[0].x = c;
4626858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
4636858538eSMatthew G. Knepley   PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL));
4646858538eSMatthew G. Knepley   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL));
4656858538eSMatthew G. Knepley   PetscFunctionReturn(0);
4666858538eSMatthew G. Knepley }
4676858538eSMatthew G. Knepley 
4686858538eSMatthew G. Knepley /*@
469*dce8aebaSBarry Smith   DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`.
4706858538eSMatthew G. Knepley 
4716858538eSMatthew G. Knepley   Collective on dm
4726858538eSMatthew G. Knepley 
4736858538eSMatthew G. Knepley   Input Parameter:
474*dce8aebaSBarry Smith . dm - the `DM`
4756858538eSMatthew G. Knepley 
4766858538eSMatthew G. Knepley   Output Parameter:
4776858538eSMatthew G. Knepley . c - global coordinate vector
4786858538eSMatthew G. Knepley 
479*dce8aebaSBarry Smith   Level: intermediate
480*dce8aebaSBarry Smith 
481*dce8aebaSBarry Smith   Notes:
482*dce8aebaSBarry Smith   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
4836858538eSMatthew G. Knepley   destroyed the array will no longer be valid.
4846858538eSMatthew G. Knepley 
4856858538eSMatthew G. Knepley   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
4866858538eSMatthew G. Knepley 
487*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
4886858538eSMatthew G. Knepley @*/
489d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c)
490d71ae5a4SJacob Faibussowitsch {
4916858538eSMatthew G. Knepley   PetscFunctionBegin;
4926858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4936858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
4946858538eSMatthew G. Knepley   if (!dm->coordinates[1].x && dm->coordinates[1].xl) {
4956858538eSMatthew G. Knepley     DM cdm = NULL;
4966858538eSMatthew G. Knepley 
49711ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
4986858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x));
4996858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates"));
5006858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
5016858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
5026858538eSMatthew G. Knepley   }
5036858538eSMatthew G. Knepley   *c = dm->coordinates[1].x;
5046858538eSMatthew G. Knepley   PetscFunctionReturn(0);
5056858538eSMatthew G. Knepley }
5066858538eSMatthew G. Knepley 
5076858538eSMatthew G. Knepley /*@
508*dce8aebaSBarry Smith   DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates
5096858538eSMatthew G. Knepley 
5106858538eSMatthew G. Knepley   Collective on dm
5116858538eSMatthew G. Knepley 
5126858538eSMatthew G. Knepley   Input Parameters:
513*dce8aebaSBarry Smith + dm - the `DM`
5146858538eSMatthew G. Knepley - c - cellwise coordinate vector
5156858538eSMatthew G. Knepley 
5166858538eSMatthew G. Knepley   Level: intermediate
5176858538eSMatthew G. Knepley 
518*dce8aebaSBarry Smith   Notes:
519*dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
520*dce8aebaSBarry Smith 
521*dce8aebaSBarry Smith   The vector c should be destroyed by the caller.
522*dce8aebaSBarry Smith 
523*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
5246858538eSMatthew G. Knepley @*/
525d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinates(DM dm, Vec c)
526d71ae5a4SJacob Faibussowitsch {
5276858538eSMatthew G. Knepley   PetscFunctionBegin;
5286858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5296858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
5306858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
5316858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
5326858538eSMatthew G. Knepley   dm->coordinates[1].x = c;
5336858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
5346858538eSMatthew G. Knepley   PetscFunctionReturn(0);
5356858538eSMatthew G. Knepley }
5366858538eSMatthew G. Knepley 
5376858538eSMatthew G. Knepley /*@
538*dce8aebaSBarry Smith   DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
5396858538eSMatthew G. Knepley 
5406858538eSMatthew G. Knepley   Collective on dm
5416858538eSMatthew G. Knepley 
5426858538eSMatthew G. Knepley   Input Parameter:
543*dce8aebaSBarry Smith . dm - the `DM`
5446858538eSMatthew G. Knepley 
5456858538eSMatthew G. Knepley   Level: advanced
5466858538eSMatthew G. Knepley 
547*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMGetCoordinatesLocalNoncollective()`
5486858538eSMatthew G. Knepley @*/
549d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm)
550d71ae5a4SJacob Faibussowitsch {
5516858538eSMatthew G. Knepley   PetscFunctionBegin;
5526858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5536858538eSMatthew G. Knepley   if (!dm->coordinates[0].xl && dm->coordinates[0].x) {
5546858538eSMatthew G. Knepley     DM cdm = NULL;
5556858538eSMatthew G. Knepley 
5566858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
5576858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
5586858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
5596858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5606858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5616858538eSMatthew G. Knepley   }
5626858538eSMatthew G. Knepley   PetscFunctionReturn(0);
5636858538eSMatthew G. Knepley }
5646858538eSMatthew G. Knepley 
5656858538eSMatthew G. Knepley /*@
566*dce8aebaSBarry Smith   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.
5676858538eSMatthew G. Knepley 
568*dce8aebaSBarry Smith   Collective on dm the first time it is called
5696858538eSMatthew G. Knepley 
5706858538eSMatthew G. Knepley   Input Parameter:
571*dce8aebaSBarry Smith . dm - the `DM`
5726858538eSMatthew G. Knepley 
5736858538eSMatthew G. Knepley   Output Parameter:
5746858538eSMatthew G. Knepley . c - coordinate vector
5756858538eSMatthew G. Knepley 
576*dce8aebaSBarry Smith   Level: intermediate
577*dce8aebaSBarry Smith 
578*dce8aebaSBarry Smith   Notes:
5796858538eSMatthew G. Knepley   This is a borrowed reference, so the user should NOT destroy this vector
5806858538eSMatthew G. Knepley 
5816858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
5826858538eSMatthew G. Knepley 
583*dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
5846858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
5856858538eSMatthew G. Knepley 
586*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
5876858538eSMatthew G. Knepley @*/
588d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
589d71ae5a4SJacob Faibussowitsch {
5906858538eSMatthew G. Knepley   PetscFunctionBegin;
5916858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5926858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
5936858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
5946858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
5956858538eSMatthew G. Knepley   PetscFunctionReturn(0);
5966858538eSMatthew G. Knepley }
5976858538eSMatthew G. Knepley 
5986858538eSMatthew G. Knepley /*@
599*dce8aebaSBarry Smith   DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.
6006858538eSMatthew G. Knepley 
6016858538eSMatthew G. Knepley   Not collective
6026858538eSMatthew G. Knepley 
6036858538eSMatthew G. Knepley   Input Parameter:
604*dce8aebaSBarry Smith . dm - the `DM`
6056858538eSMatthew G. Knepley 
6066858538eSMatthew G. Knepley   Output Parameter:
6076858538eSMatthew G. Knepley . c - coordinate vector
6086858538eSMatthew G. Knepley 
6096858538eSMatthew G. Knepley   Level: advanced
6106858538eSMatthew G. Knepley 
611*dce8aebaSBarry Smith   Note:
612*dce8aebaSBarry Smith   A previous call to  `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.
613*dce8aebaSBarry Smith 
614*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6156858538eSMatthew G. Knepley @*/
616d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
617d71ae5a4SJacob Faibussowitsch {
6186858538eSMatthew G. Knepley   PetscFunctionBegin;
6196858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6206858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
6216858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
6226858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
6236858538eSMatthew G. Knepley   PetscFunctionReturn(0);
6246858538eSMatthew G. Knepley }
6256858538eSMatthew G. Knepley 
6266858538eSMatthew G. Knepley /*@
627*dce8aebaSBarry Smith   DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.
6286858538eSMatthew G. Knepley 
6296858538eSMatthew G. Knepley   Not collective
6306858538eSMatthew G. Knepley 
6316858538eSMatthew G. Knepley   Input Parameters:
632*dce8aebaSBarry Smith + dm - the `DM`
633*dce8aebaSBarry Smith - p - the `IS` of points whose coordinates will be returned
6346858538eSMatthew G. Knepley 
6356858538eSMatthew G. Knepley   Output Parameters:
636*dce8aebaSBarry Smith + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in p, and DOFs correspond to coordinates
637*dce8aebaSBarry Smith - pCoord - the `Vec` with coordinates of points in p
6386858538eSMatthew G. Knepley 
639*dce8aebaSBarry Smith   Level: advanced
640*dce8aebaSBarry Smith 
641*dce8aebaSBarry Smith   Notes:
642*dce8aebaSBarry Smith   `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.
6436858538eSMatthew G. Knepley 
6446858538eSMatthew G. Knepley   This creates a new vector, so the user SHOULD destroy this vector
6456858538eSMatthew G. Knepley 
6466858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
6476858538eSMatthew G. Knepley 
648*dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
6496858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
6506858538eSMatthew G. Knepley 
651*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6526858538eSMatthew G. Knepley @*/
653d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
654d71ae5a4SJacob Faibussowitsch {
6556858538eSMatthew G. Knepley   DM                 cdm;
6566858538eSMatthew G. Knepley   PetscSection       cs, newcs;
6576858538eSMatthew G. Knepley   Vec                coords;
6586858538eSMatthew G. Knepley   const PetscScalar *arr;
6596858538eSMatthew G. Knepley   PetscScalar       *newarr = NULL;
6606858538eSMatthew G. Knepley   PetscInt           n;
6616858538eSMatthew G. Knepley 
6626858538eSMatthew G. Knepley   PetscFunctionBegin;
6636858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6646858538eSMatthew G. Knepley   PetscValidHeaderSpecific(p, IS_CLASSID, 2);
6656858538eSMatthew G. Knepley   if (pCoordSection) PetscValidPointer(pCoordSection, 3);
6666858538eSMatthew G. Knepley   if (pCoord) PetscValidPointer(pCoord, 4);
6676858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
6686858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, &cs));
6696858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
6706858538eSMatthew G. Knepley   PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
6716858538eSMatthew G. Knepley   PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
6726858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(coords, &arr));
6736858538eSMatthew G. Knepley   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
6746858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(coords, &arr));
6756858538eSMatthew G. Knepley   if (pCoord) {
6766858538eSMatthew G. Knepley     PetscCall(PetscSectionGetStorageSize(newcs, &n));
6776858538eSMatthew G. Knepley     /* set array in two steps to mimic PETSC_OWN_POINTER */
6786858538eSMatthew G. Knepley     PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
6796858538eSMatthew G. Knepley     PetscCall(VecReplaceArray(*pCoord, newarr));
6806858538eSMatthew G. Knepley   } else {
6816858538eSMatthew G. Knepley     PetscCall(PetscFree(newarr));
6826858538eSMatthew G. Knepley   }
6839371c9d4SSatish Balay   if (pCoordSection) {
6849371c9d4SSatish Balay     *pCoordSection = newcs;
6859371c9d4SSatish Balay   } else PetscCall(PetscSectionDestroy(&newcs));
6866858538eSMatthew G. Knepley   PetscFunctionReturn(0);
6876858538eSMatthew G. Knepley }
6886858538eSMatthew G. Knepley 
6896858538eSMatthew G. Knepley /*@
690*dce8aebaSBarry Smith   DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates
6916858538eSMatthew G. Knepley 
6926858538eSMatthew G. Knepley   Not collective
6936858538eSMatthew G. Knepley 
6946858538eSMatthew G. Knepley    Input Parameters:
695*dce8aebaSBarry Smith +  dm - the `DM`
6966858538eSMatthew G. Knepley -  c - coordinate vector
6976858538eSMatthew G. Knepley 
698*dce8aebaSBarry Smith   Level: intermediate
699*dce8aebaSBarry Smith 
7006858538eSMatthew G. Knepley   Notes:
701*dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
702*dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
7036858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
7046858538eSMatthew G. Knepley 
7056858538eSMatthew G. Knepley   The vector c should be destroyed by the caller.
7066858538eSMatthew G. Knepley 
707*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
7086858538eSMatthew G. Knepley @*/
709d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
710d71ae5a4SJacob Faibussowitsch {
7116858538eSMatthew G. Knepley   PetscFunctionBegin;
7126858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7136858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
7146858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
7156858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
7166858538eSMatthew G. Knepley   dm->coordinates[0].xl = c;
7176858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
7186858538eSMatthew G. Knepley   PetscFunctionReturn(0);
7196858538eSMatthew G. Knepley }
7206858538eSMatthew G. Knepley 
7216858538eSMatthew G. Knepley /*@
722*dce8aebaSBarry Smith   DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
7236858538eSMatthew G. Knepley 
7246858538eSMatthew G. Knepley   Collective on dm
7256858538eSMatthew G. Knepley 
7266858538eSMatthew G. Knepley   Input Parameter:
727*dce8aebaSBarry Smith . dm - the `DM`
7286858538eSMatthew G. Knepley 
7296858538eSMatthew G. Knepley   Level: advanced
7306858538eSMatthew G. Knepley 
731*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
7326858538eSMatthew G. Knepley @*/
733d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
734d71ae5a4SJacob Faibussowitsch {
7356858538eSMatthew G. Knepley   PetscFunctionBegin;
7366858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7376858538eSMatthew G. Knepley   if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
7386858538eSMatthew G. Knepley     DM cdm = NULL;
7396858538eSMatthew G. Knepley 
74011ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
7416858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
7426858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
7436858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7446858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7456858538eSMatthew G. Knepley   }
7466858538eSMatthew G. Knepley   PetscFunctionReturn(0);
7476858538eSMatthew G. Knepley }
7486858538eSMatthew G. Knepley 
7496858538eSMatthew G. Knepley /*@
750*dce8aebaSBarry Smith   DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.
7516858538eSMatthew G. Knepley 
7526858538eSMatthew G. Knepley   Collective on dm
7536858538eSMatthew G. Knepley 
7546858538eSMatthew G. Knepley   Input Parameter:
755*dce8aebaSBarry Smith . dm - the `DM`
7566858538eSMatthew G. Knepley 
7576858538eSMatthew G. Knepley   Output Parameter:
7586858538eSMatthew G. Knepley . c - coordinate vector
7596858538eSMatthew G. Knepley 
760*dce8aebaSBarry Smith   Level: intermediate
761*dce8aebaSBarry Smith 
762*dce8aebaSBarry Smith   Notes:
7636858538eSMatthew G. Knepley   This is a borrowed reference, so the user should NOT destroy this vector
7646858538eSMatthew G. Knepley 
7656858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
7666858538eSMatthew G. Knepley 
767*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
7686858538eSMatthew G. Knepley @*/
769d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
770d71ae5a4SJacob Faibussowitsch {
7716858538eSMatthew G. Knepley   PetscFunctionBegin;
7726858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7736858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
7746858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
7756858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
7766858538eSMatthew G. Knepley   PetscFunctionReturn(0);
7776858538eSMatthew G. Knepley }
7786858538eSMatthew G. Knepley 
7796858538eSMatthew G. Knepley /*@
780*dce8aebaSBarry Smith   DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.
7816858538eSMatthew G. Knepley 
7826858538eSMatthew G. Knepley   Not collective
7836858538eSMatthew G. Knepley 
7846858538eSMatthew G. Knepley   Input Parameter:
785*dce8aebaSBarry Smith . dm - the `DM`
7866858538eSMatthew G. Knepley 
7876858538eSMatthew G. Knepley   Output Parameter:
7886858538eSMatthew G. Knepley . c - cellwise coordinate vector
7896858538eSMatthew G. Knepley 
7906858538eSMatthew G. Knepley   Level: advanced
7916858538eSMatthew G. Knepley 
792*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
7936858538eSMatthew G. Knepley @*/
794d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
795d71ae5a4SJacob Faibussowitsch {
7966858538eSMatthew G. Knepley   PetscFunctionBegin;
7976858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7986858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
7996858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
8006858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
8016858538eSMatthew G. Knepley   PetscFunctionReturn(0);
8026858538eSMatthew G. Knepley }
8036858538eSMatthew G. Knepley 
8046858538eSMatthew G. Knepley /*@
805*dce8aebaSBarry Smith   DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates
8066858538eSMatthew G. Knepley 
8076858538eSMatthew G. Knepley   Not collective
8086858538eSMatthew G. Knepley 
8096858538eSMatthew G. Knepley    Input Parameters:
810*dce8aebaSBarry Smith +  dm - the `DM`
8116858538eSMatthew G. Knepley -  c - cellwise coordinate vector
8126858538eSMatthew G. Knepley 
813*dce8aebaSBarry Smith   Level: intermediate
814*dce8aebaSBarry Smith 
8156858538eSMatthew G. Knepley   Notes:
816*dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
817*dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
8186858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
8196858538eSMatthew G. Knepley 
8206858538eSMatthew G. Knepley   The vector c should be destroyed by the caller.
8216858538eSMatthew G. Knepley 
822*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
8236858538eSMatthew G. Knepley @*/
824d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
825d71ae5a4SJacob Faibussowitsch {
8266858538eSMatthew G. Knepley   PetscFunctionBegin;
8276858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8286858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
8296858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
8306858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
8316858538eSMatthew G. Knepley   dm->coordinates[1].xl = c;
8326858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
8336858538eSMatthew G. Knepley   PetscFunctionReturn(0);
8346858538eSMatthew G. Knepley }
8356858538eSMatthew G. Knepley 
836d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
837d71ae5a4SJacob Faibussowitsch {
8386858538eSMatthew G. Knepley   PetscFunctionBegin;
8396858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8406858538eSMatthew G. Knepley   PetscValidPointer(field, 2);
8416858538eSMatthew G. Knepley   if (!dm->coordinates[0].field) {
8426858538eSMatthew G. Knepley     if (dm->ops->createcoordinatefield) PetscCall((*dm->ops->createcoordinatefield)(dm, &dm->coordinates[0].field));
8436858538eSMatthew G. Knepley   }
8446858538eSMatthew G. Knepley   *field = dm->coordinates[0].field;
8456858538eSMatthew G. Knepley   PetscFunctionReturn(0);
8466858538eSMatthew G. Knepley }
8476858538eSMatthew G. Knepley 
848d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
849d71ae5a4SJacob Faibussowitsch {
8506858538eSMatthew G. Knepley   PetscFunctionBegin;
8516858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8526858538eSMatthew G. Knepley   if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2);
8536858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)field));
8546858538eSMatthew G. Knepley   PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
8556858538eSMatthew G. Knepley   dm->coordinates[0].field = field;
8566858538eSMatthew G. Knepley   PetscFunctionReturn(0);
8576858538eSMatthew G. Knepley }
8586858538eSMatthew G. Knepley 
8596858538eSMatthew G. Knepley /*@
860*dce8aebaSBarry Smith   DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.
8616858538eSMatthew G. Knepley 
8626858538eSMatthew G. Knepley   Not collective
8636858538eSMatthew G. Knepley 
8646858538eSMatthew G. Knepley   Input Parameter:
865*dce8aebaSBarry Smith . dm - the `DM`
8666858538eSMatthew G. Knepley 
8676858538eSMatthew G. Knepley   Output Parameters:
8686858538eSMatthew G. Knepley + lmin - local minimum coordinates (length coord dim, optional)
8696858538eSMatthew G. Knepley - lmax - local maximim coordinates (length coord dim, optional)
8706858538eSMatthew G. Knepley 
8716858538eSMatthew G. Knepley   Level: beginner
8726858538eSMatthew G. Knepley 
873*dce8aebaSBarry Smith   Note:
874*dce8aebaSBarry Smith   If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.
8756858538eSMatthew G. Knepley 
876*dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
8776858538eSMatthew G. Knepley @*/
878d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
879d71ae5a4SJacob Faibussowitsch {
8806858538eSMatthew G. Knepley   Vec       coords = NULL;
8816858538eSMatthew G. Knepley   PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
8826858538eSMatthew G. Knepley   PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
8836858538eSMatthew G. Knepley   PetscInt  cdim, i, j;
8846858538eSMatthew G. Knepley 
8856858538eSMatthew G. Knepley   PetscFunctionBegin;
8866858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8876858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
888bdf15c88SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
889bdf15c88SMatthew G. Knepley   if (coords) {
890bdf15c88SMatthew G. Knepley     const PetscScalar *local_coords;
891bdf15c88SMatthew G. Knepley     PetscInt           N, Ni;
892bdf15c88SMatthew G. Knepley 
893bdf15c88SMatthew G. Knepley     for (j = cdim; j < 3; ++j) {
894bdf15c88SMatthew G. Knepley       min[j] = 0;
895bdf15c88SMatthew G. Knepley       max[j] = 0;
896bdf15c88SMatthew G. Knepley     }
897bdf15c88SMatthew G. Knepley     PetscCall(VecGetArrayRead(coords, &local_coords));
898bdf15c88SMatthew G. Knepley     PetscCall(VecGetLocalSize(coords, &N));
899bdf15c88SMatthew G. Knepley     Ni = N / cdim;
900bdf15c88SMatthew G. Knepley     for (i = 0; i < Ni; ++i) {
901bdf15c88SMatthew G. Knepley       for (j = 0; j < cdim; ++j) {
902bdf15c88SMatthew G. Knepley         min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
903bdf15c88SMatthew G. Knepley         max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
904bdf15c88SMatthew G. Knepley       }
905bdf15c88SMatthew G. Knepley     }
906bdf15c88SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(coords, &local_coords));
907bdf15c88SMatthew G. Knepley     PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
9086858538eSMatthew G. Knepley     if (coords) {
9096858538eSMatthew G. Knepley       PetscCall(VecGetArrayRead(coords, &local_coords));
9106858538eSMatthew G. Knepley       PetscCall(VecGetLocalSize(coords, &N));
9116858538eSMatthew G. Knepley       Ni = N / cdim;
9126858538eSMatthew G. Knepley       for (i = 0; i < Ni; ++i) {
913bdf15c88SMatthew G. Knepley         for (j = 0; j < cdim; ++j) {
914bdf15c88SMatthew G. Knepley           min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
915bdf15c88SMatthew G. Knepley           max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
9166858538eSMatthew G. Knepley         }
9176858538eSMatthew G. Knepley       }
9186858538eSMatthew G. Knepley       PetscCall(VecRestoreArrayRead(coords, &local_coords));
919bdf15c88SMatthew G. Knepley     }
9206858538eSMatthew G. Knepley   } else {
9216858538eSMatthew G. Knepley     PetscBool isda;
9226858538eSMatthew G. Knepley 
9236858538eSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
9246858538eSMatthew G. Knepley     if (isda) PetscCall(DMGetLocalBoundingIndices_DMDA(dm, min, max));
9256858538eSMatthew G. Knepley   }
9266858538eSMatthew G. Knepley   if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
9276858538eSMatthew G. Knepley   if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
9286858538eSMatthew G. Knepley   PetscFunctionReturn(0);
9296858538eSMatthew G. Knepley }
9306858538eSMatthew G. Knepley 
9316858538eSMatthew G. Knepley /*@
932*dce8aebaSBarry Smith   DMGetBoundingBox - Returns the global bounding box for the `DM`.
9336858538eSMatthew G. Knepley 
9346858538eSMatthew G. Knepley   Collective
9356858538eSMatthew G. Knepley 
9366858538eSMatthew G. Knepley   Input Parameter:
937*dce8aebaSBarry Smith . dm - the `DM`
9386858538eSMatthew G. Knepley 
9396858538eSMatthew G. Knepley   Output Parameters:
9406858538eSMatthew G. Knepley + gmin - global minimum coordinates (length coord dim, optional)
9416858538eSMatthew G. Knepley - gmax - global maximim coordinates (length coord dim, optional)
9426858538eSMatthew G. Knepley 
9436858538eSMatthew G. Knepley   Level: beginner
9446858538eSMatthew G. Knepley 
945*dce8aebaSBarry Smith .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
9466858538eSMatthew G. Knepley @*/
947d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
948d71ae5a4SJacob Faibussowitsch {
9496858538eSMatthew G. Knepley   PetscReal   lmin[3], lmax[3];
9506858538eSMatthew G. Knepley   PetscInt    cdim;
9516858538eSMatthew G. Knepley   PetscMPIInt count;
9526858538eSMatthew G. Knepley 
9536858538eSMatthew G. Knepley   PetscFunctionBegin;
9546858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9556858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
9566858538eSMatthew G. Knepley   PetscCall(PetscMPIIntCast(cdim, &count));
9576858538eSMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
9586858538eSMatthew G. Knepley   if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
9596858538eSMatthew G. Knepley   if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
9606858538eSMatthew G. Knepley   PetscFunctionReturn(0);
9616858538eSMatthew G. Knepley }
9626858538eSMatthew G. Knepley 
9636858538eSMatthew G. Knepley /*@
9646858538eSMatthew G. Knepley   DMProjectCoordinates - Project coordinates to a different space
9656858538eSMatthew G. Knepley 
9666858538eSMatthew G. Knepley   Input Parameters:
967*dce8aebaSBarry Smith + dm      - The `DM` object
9686858538eSMatthew G. Knepley - disc    - The new coordinate discretization or NULL to ensure a coordinate discretization exists
9696858538eSMatthew G. Knepley 
9706858538eSMatthew G. Knepley   Level: intermediate
9716858538eSMatthew G. Knepley 
972*dce8aebaSBarry Smith   Notes:
973*dce8aebaSBarry Smith   A `PetscFE` defines an approximation space using a `PetscSpace`, which represents the basis functions, and a `PetscDualSpace`, which defines the interpolation operation
974*dce8aebaSBarry Smith   in the space.
975*dce8aebaSBarry Smith 
976*dce8aebaSBarry Smith   This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
977*dce8aebaSBarry Smith   The coordinate projection is done on the continuous coordinates, and if possible, the discontinuous coordinates are also updated.
978*dce8aebaSBarry Smith 
979*dce8aebaSBarry Smith   Developer Note:
980*dce8aebaSBarry Smith   With more effort, we could directly project the discontinuous coordinates also.
981*dce8aebaSBarry Smith 
982*dce8aebaSBarry Smith .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
9836858538eSMatthew G. Knepley @*/
984d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectCoordinates(DM dm, PetscFE disc)
985d71ae5a4SJacob Faibussowitsch {
9866858538eSMatthew G. Knepley   PetscFE      discOld;
9876858538eSMatthew G. Knepley   PetscClassId classid;
9886858538eSMatthew G. Knepley   DM           cdmOld, cdmNew;
9896858538eSMatthew G. Knepley   Vec          coordsOld, coordsNew;
9906858538eSMatthew G. Knepley   Mat          matInterp;
9916858538eSMatthew G. Knepley   PetscBool    same_space = PETSC_TRUE;
9926858538eSMatthew G. Knepley 
9936858538eSMatthew G. Knepley   PetscFunctionBegin;
9946858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9956858538eSMatthew G. Knepley   if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2);
9966858538eSMatthew G. Knepley 
9976858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdmOld));
9986858538eSMatthew G. Knepley   /* Check current discretization is compatible */
9996858538eSMatthew G. Knepley   PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
10006858538eSMatthew G. Knepley   PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
10016858538eSMatthew G. Knepley   if (classid != PETSCFE_CLASSID) {
10026858538eSMatthew G. Knepley     if (classid == PETSC_CONTAINER_CLASSID) {
10036858538eSMatthew G. Knepley       PetscFE        feLinear;
10046858538eSMatthew G. Knepley       DMPolytopeType ct;
10056858538eSMatthew G. Knepley       PetscInt       dim, dE, cStart, cEnd;
10066858538eSMatthew G. Knepley       PetscBool      simplex;
10076858538eSMatthew G. Knepley 
10086858538eSMatthew G. Knepley       /* Assume linear vertex coordinates */
10096858538eSMatthew G. Knepley       PetscCall(DMGetDimension(dm, &dim));
10106858538eSMatthew G. Knepley       PetscCall(DMGetCoordinateDim(dm, &dE));
10116858538eSMatthew G. Knepley       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
10126858538eSMatthew G. Knepley       if (cEnd > cStart) {
10136858538eSMatthew G. Knepley         PetscCall(DMPlexGetCellType(dm, cStart, &ct));
10146858538eSMatthew G. Knepley         switch (ct) {
10156858538eSMatthew G. Knepley         case DM_POLYTOPE_TRI_PRISM:
1016d71ae5a4SJacob Faibussowitsch         case DM_POLYTOPE_TRI_PRISM_TENSOR:
1017d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot autoamtically create coordinate space for prisms");
1018d71ae5a4SJacob Faibussowitsch         default:
1019d71ae5a4SJacob Faibussowitsch           break;
10206858538eSMatthew G. Knepley         }
10216858538eSMatthew G. Knepley       }
10226858538eSMatthew G. Knepley       PetscCall(DMPlexIsSimplex(dm, &simplex));
10236858538eSMatthew G. Knepley       PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dE, simplex, 1, -1, &feLinear));
10246858538eSMatthew G. Knepley       PetscCall(DMSetField(cdmOld, 0, NULL, (PetscObject)feLinear));
10256858538eSMatthew G. Knepley       PetscCall(PetscFEDestroy(&feLinear));
10266858538eSMatthew G. Knepley       PetscCall(DMCreateDS(cdmOld));
10276858538eSMatthew G. Knepley       PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
10286858538eSMatthew G. Knepley     } else {
10296858538eSMatthew G. Knepley       const char *discname;
10306858538eSMatthew G. Knepley 
10316858538eSMatthew G. Knepley       PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
10326858538eSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
10336858538eSMatthew G. Knepley     }
10346858538eSMatthew G. Knepley   }
10356858538eSMatthew G. Knepley   if (!disc) PetscFunctionReturn(0);
10366858538eSMatthew G. Knepley   { // Check if the new space is the same as the old modulo quadrature
10376858538eSMatthew G. Knepley     PetscDualSpace dsOld, ds;
10386858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
10396858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(disc, &ds));
10406858538eSMatthew G. Knepley     PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
10416858538eSMatthew G. Knepley   }
10426858538eSMatthew G. Knepley   /* Make a fresh clone of the coordinate DM */
10436858538eSMatthew G. Knepley   PetscCall(DMClone(cdmOld, &cdmNew));
10446858538eSMatthew G. Knepley   PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
10456858538eSMatthew G. Knepley   PetscCall(DMCreateDS(cdmNew));
10466858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coordsOld));
10476858538eSMatthew G. Knepley   if (same_space) {
10486858538eSMatthew G. Knepley     PetscCall(PetscObjectReference((PetscObject)coordsOld));
10496858538eSMatthew G. Knepley     coordsNew = coordsOld;
10506858538eSMatthew G. Knepley   } else { // Project the coordinate vector from old to new space
10516858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
10526858538eSMatthew G. Knepley     PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &matInterp, NULL));
10536858538eSMatthew G. Knepley     PetscCall(MatInterpolate(matInterp, coordsOld, coordsNew));
10546858538eSMatthew G. Knepley     PetscCall(MatDestroy(&matInterp));
10556858538eSMatthew G. Knepley   }
10566858538eSMatthew G. Knepley   /* Set new coordinate structures */
10576858538eSMatthew G. Knepley   PetscCall(DMSetCoordinateField(dm, NULL));
10586858538eSMatthew G. Knepley   PetscCall(DMSetCoordinateDM(dm, cdmNew));
10596858538eSMatthew G. Knepley   PetscCall(DMSetCoordinates(dm, coordsNew));
10606858538eSMatthew G. Knepley   PetscCall(VecDestroy(&coordsNew));
10616858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cdmNew));
10626858538eSMatthew G. Knepley   PetscFunctionReturn(0);
10636858538eSMatthew G. Knepley }
10646858538eSMatthew G. Knepley 
10656858538eSMatthew G. Knepley /*@
1066*dce8aebaSBarry Smith   DMLocatePoints - Locate the points in v in the mesh and return a `PetscSF` of the containing cells
10676858538eSMatthew G. Knepley 
10686858538eSMatthew G. Knepley   Collective on v (see explanation below)
10696858538eSMatthew G. Knepley 
10706858538eSMatthew G. Knepley   Input Parameters:
1071*dce8aebaSBarry Smith + dm - The `DM`
1072*dce8aebaSBarry Smith - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`
10736858538eSMatthew G. Knepley 
10746858538eSMatthew G. Knepley   Input/Output Parameters:
1075*dce8aebaSBarry Smith + v - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
1076*dce8aebaSBarry Smith - cellSF - Points to either NULL, or a `PetscSF` with guesses for which cells contain each point;
1077*dce8aebaSBarry Smith            on output, the `PetscSF` containing the ranks and local indices of the containing points
10786858538eSMatthew G. Knepley 
10796858538eSMatthew G. Knepley   Level: developer
10806858538eSMatthew G. Knepley 
10816858538eSMatthew G. Knepley   Notes:
1082*dce8aebaSBarry Smith   To do a search of the local cells of the mesh, v should have `PETSC_COMM_SELF` as its communicator.
10836858538eSMatthew G. Knepley   To do a search of all the cells in the distributed mesh, v should have the same communicator as dm.
10846858538eSMatthew G. Knepley 
1085d8206211SMatthew G. Knepley   Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.
1086d8206211SMatthew G. Knepley 
1087*dce8aebaSBarry Smith   If *cellSF is NULL on input, a `PetscSF` will be created.
1088*dce8aebaSBarry Smith   If *cellSF is not NULL on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.
10896858538eSMatthew G. Knepley 
10906858538eSMatthew G. Knepley   An array that maps each point to its containing cell can be obtained with
1091*dce8aebaSBarry Smith .vb
1092*dce8aebaSBarry Smith     const PetscSFNode *cells;
1093*dce8aebaSBarry Smith     PetscInt           nFound;
1094*dce8aebaSBarry Smith     const PetscInt    *found;
10956858538eSMatthew G. Knepley 
1096*dce8aebaSBarry Smith     PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1097*dce8aebaSBarry Smith .ve
10986858538eSMatthew G. Knepley 
10996858538eSMatthew G. Knepley   Where cells[i].rank is the rank of the cell containing point found[i] (or i if found == NULL), and cells[i].index is
11006858538eSMatthew G. Knepley   the index of the cell in its rank's local numbering.
11016858538eSMatthew G. Knepley 
1102*dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
11036858538eSMatthew G. Knepley @*/
1104d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1105d71ae5a4SJacob Faibussowitsch {
11066858538eSMatthew G. Knepley   PetscFunctionBegin;
11076858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11086858538eSMatthew G. Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 2);
11096858538eSMatthew G. Knepley   PetscValidPointer(cellSF, 4);
11106858538eSMatthew G. Knepley   if (*cellSF) {
11116858538eSMatthew G. Knepley     PetscMPIInt result;
11126858538eSMatthew G. Knepley 
11136858538eSMatthew G. Knepley     PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4);
11146858538eSMatthew G. Knepley     PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
11156858538eSMatthew 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");
11166858538eSMatthew G. Knepley   } else {
11176858538eSMatthew G. Knepley     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
11186858538eSMatthew G. Knepley   }
11196858538eSMatthew G. Knepley   PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1120dbbe0bcdSBarry Smith   PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
11216858538eSMatthew G. Knepley   PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
11226858538eSMatthew G. Knepley   PetscFunctionReturn(0);
11236858538eSMatthew G. Knepley }
1124