xref: /petsc/src/dm/interface/dmcoordinates.c (revision 20f4b53cbb5e9bd9ef12b76a8697d60d197cda17)
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   }
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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   }
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
586858538eSMatthew G. Knepley }
596858538eSMatthew G. Knepley 
606858538eSMatthew G. Knepley /*@
61dce8aebaSBarry Smith   DMGetCoordinateDM - Gets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
626858538eSMatthew G. Knepley 
63*20f4b53cSBarry Smith   Collective
646858538eSMatthew G. Knepley 
656858538eSMatthew G. Knepley   Input Parameter:
66dce8aebaSBarry Smith . dm - the `DM`
676858538eSMatthew G. Knepley 
686858538eSMatthew G. Knepley   Output Parameter:
69dce8aebaSBarry Smith . cdm - coordinate `DM`
706858538eSMatthew G. Knepley 
716858538eSMatthew G. Knepley   Level: intermediate
726858538eSMatthew G. Knepley 
73dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGSetCellCoordinateDM()`,
74dce8aebaSBarry 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;
923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
936858538eSMatthew G. Knepley }
946858538eSMatthew G. Knepley 
956858538eSMatthew G. Knepley /*@
96dce8aebaSBarry Smith   DMSetCoordinateDM - Sets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
976858538eSMatthew G. Knepley 
98*20f4b53cSBarry Smith   Logically Collective
996858538eSMatthew G. Knepley 
1006858538eSMatthew G. Knepley   Input Parameters:
101dce8aebaSBarry Smith + dm - the `DM`
102dce8aebaSBarry Smith - cdm - coordinate `DM`
1036858538eSMatthew G. Knepley 
1046858538eSMatthew G. Knepley   Level: intermediate
1056858538eSMatthew G. Knepley 
106dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMGetCellCoordinateDM()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`,
107dce8aebaSBarry 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);
11390612307SJed Brown   if (cdm) 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;
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1186858538eSMatthew G. Knepley }
1196858538eSMatthew G. Knepley 
1206858538eSMatthew G. Knepley /*@
121dce8aebaSBarry Smith   DMGetCellCoordinateDM - Gets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1226858538eSMatthew G. Knepley 
123*20f4b53cSBarry Smith   Collective
1246858538eSMatthew G. Knepley 
1256858538eSMatthew G. Knepley   Input Parameter:
126dce8aebaSBarry Smith . dm - the `DM`
1276858538eSMatthew G. Knepley 
1286858538eSMatthew G. Knepley   Output Parameter:
129dce8aebaSBarry Smith . cdm - cellwise coordinate `DM`, or NULL if they are not defined
1306858538eSMatthew G. Knepley 
1316858538eSMatthew G. Knepley   Level: intermediate
1326858538eSMatthew G. Knepley 
133dce8aebaSBarry Smith   Note:
134dce8aebaSBarry Smith   Call `DMLocalizeCoordinates()` to automatically create cellwise coordinates for periodic geometries.
135dce8aebaSBarry Smith 
136dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
137dce8aebaSBarry 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;
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1466858538eSMatthew G. Knepley }
1476858538eSMatthew G. Knepley 
1486858538eSMatthew G. Knepley /*@
149dce8aebaSBarry Smith   DMSetCellCoordinateDM - Sets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1506858538eSMatthew G. Knepley 
151*20f4b53cSBarry Smith   Logically Collective
1526858538eSMatthew G. Knepley 
1536858538eSMatthew G. Knepley   Input Parameters:
154dce8aebaSBarry Smith + dm - the `DM`
155dce8aebaSBarry Smith - cdm - cellwise coordinate `DM`
1566858538eSMatthew G. Knepley 
1576858538eSMatthew G. Knepley   Level: intermediate
1586858538eSMatthew G. Knepley 
159dce8aebaSBarry Smith   Note:
16035cb6cd3SPierre Jolivet   As opposed to `DMSetCoordinateDM()` these coordinates are useful for discontinuous Galerkin methods since they support coordinate fields that are discontinuous at cell boundaries.
161dce8aebaSBarry Smith 
162dce8aebaSBarry Smith .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
163dce8aebaSBarry 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;
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1806858538eSMatthew G. Knepley }
1816858538eSMatthew G. Knepley 
1826858538eSMatthew G. Knepley /*@
183dce8aebaSBarry 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:
188dce8aebaSBarry 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 
195dce8aebaSBarry 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;
2043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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:
213dce8aebaSBarry Smith + dm  - The `DM` object
2146858538eSMatthew G. Knepley - dim - The embedding dimension
2156858538eSMatthew G. Knepley 
2166858538eSMatthew G. Knepley   Level: intermediate
2176858538eSMatthew G. Knepley 
218dce8aebaSBarry 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   }
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
241*20f4b53cSBarry Smith   Collective
2426858538eSMatthew G. Knepley 
2436858538eSMatthew G. Knepley   Input Parameter:
244dce8aebaSBarry Smith . dm - The `DM` object
2456858538eSMatthew G. Knepley 
2466858538eSMatthew G. Knepley   Output Parameter:
247dce8aebaSBarry Smith . section - The `PetscSection` object
2486858538eSMatthew G. Knepley 
2496858538eSMatthew G. Knepley   Level: intermediate
2506858538eSMatthew G. Knepley 
251dce8aebaSBarry Smith   Note:
252dce8aebaSBarry Smith   This just retrieves the local section from the coordinate `DM`. In other words,
253dce8aebaSBarry Smith .vb
254dce8aebaSBarry Smith   DMGetCoordinateDM(dm, &cdm);
255dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
256dce8aebaSBarry Smith .ve
257dce8aebaSBarry 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));
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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:
278dce8aebaSBarry Smith + dm      - The `DM` object
279dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
280dce8aebaSBarry Smith - section - The `PetscSection` object
2816858538eSMatthew G. Knepley 
2826858538eSMatthew G. Knepley   Level: intermediate
2836858538eSMatthew G. Knepley 
284dce8aebaSBarry 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   }
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
318*20f4b53cSBarry Smith   Collective
3196858538eSMatthew G. Knepley 
3206858538eSMatthew G. Knepley   Input Parameter:
321dce8aebaSBarry Smith . dm - The `DM` object
3226858538eSMatthew G. Knepley 
3236858538eSMatthew G. Knepley   Output Parameter:
324*20f4b53cSBarry Smith . section - The `PetscSection` object, or `NULL` if no cellwise coordinates are defined
3256858538eSMatthew G. Knepley 
3266858538eSMatthew G. Knepley   Level: intermediate
3276858538eSMatthew G. Knepley 
328dce8aebaSBarry Smith   Note:
329dce8aebaSBarry Smith   This just retrieves the local section from the cell coordinate `DM`. In other words,
330dce8aebaSBarry Smith .vb
331dce8aebaSBarry Smith   DMGetCellCoordinateDM(dm, &cdm);
332dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
333dce8aebaSBarry Smith .ve
334dce8aebaSBarry Smith 
335dce8aebaSBarry 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));
3473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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:
356dce8aebaSBarry Smith + dm      - The `DM` object
357dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
358dce8aebaSBarry Smith - section - The `PetscSection` object for a cellwise layout
3596858538eSMatthew G. Knepley 
3606858538eSMatthew G. Knepley   Level: intermediate
3616858538eSMatthew G. Knepley 
362dce8aebaSBarry 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   }
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3926858538eSMatthew G. Knepley }
3936858538eSMatthew G. Knepley 
3946858538eSMatthew G. Knepley /*@
395dce8aebaSBarry Smith   DMGetCoordinates - Gets a global vector with the coordinates associated with the `DM`.
3966858538eSMatthew G. Knepley 
397*20f4b53cSBarry Smith   Collective
3986858538eSMatthew G. Knepley 
3996858538eSMatthew G. Knepley   Input Parameter:
400dce8aebaSBarry Smith . dm - the `DM`
4016858538eSMatthew G. Knepley 
4026858538eSMatthew G. Knepley   Output Parameter:
4036858538eSMatthew G. Knepley . c - global coordinate vector
4046858538eSMatthew G. Knepley 
405dce8aebaSBarry Smith   Level: intermediate
406dce8aebaSBarry Smith 
407dce8aebaSBarry Smith   Notes:
408dce8aebaSBarry 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 
413dce8aebaSBarry 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 
416dce8aebaSBarry 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;
4333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4346858538eSMatthew G. Knepley }
4356858538eSMatthew G. Knepley 
4366858538eSMatthew G. Knepley /*@
437dce8aebaSBarry Smith   DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates
4386858538eSMatthew G. Knepley 
439*20f4b53cSBarry Smith   Collective
4406858538eSMatthew G. Knepley 
4416858538eSMatthew G. Knepley   Input Parameters:
442dce8aebaSBarry Smith + dm - the `DM`
4436858538eSMatthew G. Knepley - c - coordinate vector
4446858538eSMatthew G. Knepley 
4456858538eSMatthew G. Knepley   Level: intermediate
4466858538eSMatthew G. Knepley 
447dce8aebaSBarry Smith   Notes:
448dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
449dce8aebaSBarry Smith 
450*20f4b53cSBarry Smith   The vector `c` can be destroyed after the call
451dce8aebaSBarry Smith 
452dce8aebaSBarry 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));
4653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4666858538eSMatthew G. Knepley }
4676858538eSMatthew G. Knepley 
4686858538eSMatthew G. Knepley /*@
469dce8aebaSBarry Smith   DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`.
4706858538eSMatthew G. Knepley 
471*20f4b53cSBarry Smith   Collective
4726858538eSMatthew G. Knepley 
4736858538eSMatthew G. Knepley   Input Parameter:
474dce8aebaSBarry Smith . dm - the `DM`
4756858538eSMatthew G. Knepley 
4766858538eSMatthew G. Knepley   Output Parameter:
4776858538eSMatthew G. Knepley . c - global coordinate vector
4786858538eSMatthew G. Knepley 
479dce8aebaSBarry Smith   Level: intermediate
480dce8aebaSBarry Smith 
481dce8aebaSBarry Smith   Notes:
482dce8aebaSBarry 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 
487dce8aebaSBarry 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;
5043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5056858538eSMatthew G. Knepley }
5066858538eSMatthew G. Knepley 
5076858538eSMatthew G. Knepley /*@
508dce8aebaSBarry Smith   DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates
5096858538eSMatthew G. Knepley 
510*20f4b53cSBarry Smith   Collective
5116858538eSMatthew G. Knepley 
5126858538eSMatthew G. Knepley   Input Parameters:
513dce8aebaSBarry Smith + dm - the `DM`
5146858538eSMatthew G. Knepley - c - cellwise coordinate vector
5156858538eSMatthew G. Knepley 
5166858538eSMatthew G. Knepley   Level: intermediate
5176858538eSMatthew G. Knepley 
518dce8aebaSBarry Smith   Notes:
519dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
520dce8aebaSBarry Smith 
521*20f4b53cSBarry Smith   The vector `c` should be destroyed by the caller.
522dce8aebaSBarry Smith 
523dce8aebaSBarry 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));
5343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5356858538eSMatthew G. Knepley }
5366858538eSMatthew G. Knepley 
5376858538eSMatthew G. Knepley /*@
538dce8aebaSBarry Smith   DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
5396858538eSMatthew G. Knepley 
540*20f4b53cSBarry Smith   Collective
5416858538eSMatthew G. Knepley 
5426858538eSMatthew G. Knepley   Input Parameter:
543dce8aebaSBarry Smith . dm - the `DM`
5446858538eSMatthew G. Knepley 
5456858538eSMatthew G. Knepley   Level: advanced
5466858538eSMatthew G. Knepley 
547dce8aebaSBarry 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;
5559d92f5ffSMatthew G. Knepley     PetscInt bs;
5566858538eSMatthew G. Knepley 
5576858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
5586858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
5599d92f5ffSMatthew G. Knepley     // If the size of the vector is 0, it will not get the right block size
5609d92f5ffSMatthew G. Knepley     PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs));
5619d92f5ffSMatthew G. Knepley     PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs));
5626858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
5636858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5646858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5656858538eSMatthew G. Knepley   }
5663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5676858538eSMatthew G. Knepley }
5686858538eSMatthew G. Knepley 
5696858538eSMatthew G. Knepley /*@
570dce8aebaSBarry Smith   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.
5716858538eSMatthew G. Knepley 
572*20f4b53cSBarry Smith   Collective the first time it is called
5736858538eSMatthew G. Knepley 
5746858538eSMatthew G. Knepley   Input Parameter:
575dce8aebaSBarry Smith . dm - the `DM`
5766858538eSMatthew G. Knepley 
5776858538eSMatthew G. Knepley   Output Parameter:
5786858538eSMatthew G. Knepley . c - coordinate vector
5796858538eSMatthew G. Knepley 
580dce8aebaSBarry Smith   Level: intermediate
581dce8aebaSBarry Smith 
582dce8aebaSBarry Smith   Notes:
5836858538eSMatthew G. Knepley   This is a borrowed reference, so the user should NOT destroy this vector
5846858538eSMatthew G. Knepley 
5856858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
5866858538eSMatthew G. Knepley 
587dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
5886858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
5896858538eSMatthew G. Knepley 
590dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
5916858538eSMatthew G. Knepley @*/
592d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
593d71ae5a4SJacob Faibussowitsch {
5946858538eSMatthew G. Knepley   PetscFunctionBegin;
5956858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5966858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
5976858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
5986858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
5993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6006858538eSMatthew G. Knepley }
6016858538eSMatthew G. Knepley 
6026858538eSMatthew G. Knepley /*@
603dce8aebaSBarry Smith   DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.
6046858538eSMatthew G. Knepley 
605*20f4b53cSBarry Smith   Not Collective
6066858538eSMatthew G. Knepley 
6076858538eSMatthew G. Knepley   Input Parameter:
608dce8aebaSBarry Smith . dm - the `DM`
6096858538eSMatthew G. Knepley 
6106858538eSMatthew G. Knepley   Output Parameter:
6116858538eSMatthew G. Knepley . c - coordinate vector
6126858538eSMatthew G. Knepley 
6136858538eSMatthew G. Knepley   Level: advanced
6146858538eSMatthew G. Knepley 
615dce8aebaSBarry Smith   Note:
616dce8aebaSBarry Smith   A previous call to  `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.
617dce8aebaSBarry Smith 
618dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6196858538eSMatthew G. Knepley @*/
620d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
621d71ae5a4SJacob Faibussowitsch {
6226858538eSMatthew G. Knepley   PetscFunctionBegin;
6236858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6246858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
6256858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
6266858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6286858538eSMatthew G. Knepley }
6296858538eSMatthew G. Knepley 
6306858538eSMatthew G. Knepley /*@
631dce8aebaSBarry Smith   DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.
6326858538eSMatthew G. Knepley 
633*20f4b53cSBarry Smith   Not Collective
6346858538eSMatthew G. Knepley 
6356858538eSMatthew G. Knepley   Input Parameters:
636dce8aebaSBarry Smith + dm - the `DM`
637dce8aebaSBarry Smith - p - the `IS` of points whose coordinates will be returned
6386858538eSMatthew G. Knepley 
6396858538eSMatthew G. Knepley   Output Parameters:
640dce8aebaSBarry Smith + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in p, and DOFs correspond to coordinates
641dce8aebaSBarry Smith - pCoord - the `Vec` with coordinates of points in p
6426858538eSMatthew G. Knepley 
643dce8aebaSBarry Smith   Level: advanced
644dce8aebaSBarry Smith 
645dce8aebaSBarry Smith   Notes:
646dce8aebaSBarry Smith   `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.
6476858538eSMatthew G. Knepley 
6486858538eSMatthew G. Knepley   This creates a new vector, so the user SHOULD destroy this vector
6496858538eSMatthew G. Knepley 
6506858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
6516858538eSMatthew G. Knepley 
652dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
6536858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
6546858538eSMatthew G. Knepley 
655dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6566858538eSMatthew G. Knepley @*/
657d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
658d71ae5a4SJacob Faibussowitsch {
6596858538eSMatthew G. Knepley   DM                 cdm;
6606858538eSMatthew G. Knepley   PetscSection       cs, newcs;
6616858538eSMatthew G. Knepley   Vec                coords;
6626858538eSMatthew G. Knepley   const PetscScalar *arr;
6636858538eSMatthew G. Knepley   PetscScalar       *newarr = NULL;
6646858538eSMatthew G. Knepley   PetscInt           n;
6656858538eSMatthew G. Knepley 
6666858538eSMatthew G. Knepley   PetscFunctionBegin;
6676858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6686858538eSMatthew G. Knepley   PetscValidHeaderSpecific(p, IS_CLASSID, 2);
6696858538eSMatthew G. Knepley   if (pCoordSection) PetscValidPointer(pCoordSection, 3);
6706858538eSMatthew G. Knepley   if (pCoord) PetscValidPointer(pCoord, 4);
6716858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
6726858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, &cs));
6736858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
6746858538eSMatthew G. Knepley   PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
6756858538eSMatthew G. Knepley   PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
6766858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(coords, &arr));
6776858538eSMatthew G. Knepley   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
6786858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(coords, &arr));
6796858538eSMatthew G. Knepley   if (pCoord) {
6806858538eSMatthew G. Knepley     PetscCall(PetscSectionGetStorageSize(newcs, &n));
6816858538eSMatthew G. Knepley     /* set array in two steps to mimic PETSC_OWN_POINTER */
6826858538eSMatthew G. Knepley     PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
6836858538eSMatthew G. Knepley     PetscCall(VecReplaceArray(*pCoord, newarr));
6846858538eSMatthew G. Knepley   } else {
6856858538eSMatthew G. Knepley     PetscCall(PetscFree(newarr));
6866858538eSMatthew G. Knepley   }
6879371c9d4SSatish Balay   if (pCoordSection) {
6889371c9d4SSatish Balay     *pCoordSection = newcs;
6899371c9d4SSatish Balay   } else PetscCall(PetscSectionDestroy(&newcs));
6903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6916858538eSMatthew G. Knepley }
6926858538eSMatthew G. Knepley 
6936858538eSMatthew G. Knepley /*@
694dce8aebaSBarry Smith   DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates
6956858538eSMatthew G. Knepley 
696*20f4b53cSBarry Smith   Not Collective
6976858538eSMatthew G. Knepley 
6986858538eSMatthew G. Knepley    Input Parameters:
699dce8aebaSBarry Smith +  dm - the `DM`
7006858538eSMatthew G. Knepley -  c - coordinate vector
7016858538eSMatthew G. Knepley 
702dce8aebaSBarry Smith   Level: intermediate
703dce8aebaSBarry Smith 
7046858538eSMatthew G. Knepley   Notes:
705dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
706dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
7076858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
7086858538eSMatthew G. Knepley 
7096858538eSMatthew G. Knepley   The vector c should be destroyed by the caller.
7106858538eSMatthew G. Knepley 
711dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
7126858538eSMatthew G. Knepley @*/
713d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
714d71ae5a4SJacob Faibussowitsch {
7156858538eSMatthew G. Knepley   PetscFunctionBegin;
7166858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7176858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
7186858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
7196858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
7206858538eSMatthew G. Knepley   dm->coordinates[0].xl = c;
7216858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
7223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7236858538eSMatthew G. Knepley }
7246858538eSMatthew G. Knepley 
7256858538eSMatthew G. Knepley /*@
726dce8aebaSBarry Smith   DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
7276858538eSMatthew G. Knepley 
728*20f4b53cSBarry Smith   Collective
7296858538eSMatthew G. Knepley 
7306858538eSMatthew G. Knepley   Input Parameter:
731dce8aebaSBarry Smith . dm - the `DM`
7326858538eSMatthew G. Knepley 
7336858538eSMatthew G. Knepley   Level: advanced
7346858538eSMatthew G. Knepley 
735dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
7366858538eSMatthew G. Knepley @*/
737d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
738d71ae5a4SJacob Faibussowitsch {
7396858538eSMatthew G. Knepley   PetscFunctionBegin;
7406858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7416858538eSMatthew G. Knepley   if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
7426858538eSMatthew G. Knepley     DM cdm = NULL;
7436858538eSMatthew G. Knepley 
74411ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
7456858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
7466858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
7476858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7486858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7496858538eSMatthew G. Knepley   }
7503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7516858538eSMatthew G. Knepley }
7526858538eSMatthew G. Knepley 
7536858538eSMatthew G. Knepley /*@
754dce8aebaSBarry Smith   DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.
7556858538eSMatthew G. Knepley 
756*20f4b53cSBarry Smith   Collective
7576858538eSMatthew G. Knepley 
7586858538eSMatthew G. Knepley   Input Parameter:
759dce8aebaSBarry Smith . dm - the `DM`
7606858538eSMatthew G. Knepley 
7616858538eSMatthew G. Knepley   Output Parameter:
7626858538eSMatthew G. Knepley . c - coordinate vector
7636858538eSMatthew G. Knepley 
764dce8aebaSBarry Smith   Level: intermediate
765dce8aebaSBarry Smith 
766dce8aebaSBarry Smith   Notes:
7676858538eSMatthew G. Knepley   This is a borrowed reference, so the user should NOT destroy this vector
7686858538eSMatthew G. Knepley 
7696858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
7706858538eSMatthew G. Knepley 
771dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
7726858538eSMatthew G. Knepley @*/
773d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
774d71ae5a4SJacob Faibussowitsch {
7756858538eSMatthew G. Knepley   PetscFunctionBegin;
7766858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7776858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
7786858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
7796858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
7803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7816858538eSMatthew G. Knepley }
7826858538eSMatthew G. Knepley 
7836858538eSMatthew G. Knepley /*@
784dce8aebaSBarry Smith   DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.
7856858538eSMatthew G. Knepley 
786*20f4b53cSBarry Smith   Not Collective
7876858538eSMatthew G. Knepley 
7886858538eSMatthew G. Knepley   Input Parameter:
789dce8aebaSBarry Smith . dm - the `DM`
7906858538eSMatthew G. Knepley 
7916858538eSMatthew G. Knepley   Output Parameter:
7926858538eSMatthew G. Knepley . c - cellwise coordinate vector
7936858538eSMatthew G. Knepley 
7946858538eSMatthew G. Knepley   Level: advanced
7956858538eSMatthew G. Knepley 
796dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
7976858538eSMatthew G. Knepley @*/
798d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
799d71ae5a4SJacob Faibussowitsch {
8006858538eSMatthew G. Knepley   PetscFunctionBegin;
8016858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8026858538eSMatthew G. Knepley   PetscValidPointer(c, 2);
8036858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
8046858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
8053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8066858538eSMatthew G. Knepley }
8076858538eSMatthew G. Knepley 
8086858538eSMatthew G. Knepley /*@
809dce8aebaSBarry Smith   DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates
8106858538eSMatthew G. Knepley 
811*20f4b53cSBarry Smith   Not Collective
8126858538eSMatthew G. Knepley 
8136858538eSMatthew G. Knepley    Input Parameters:
814dce8aebaSBarry Smith +  dm - the `DM`
8156858538eSMatthew G. Knepley -  c - cellwise coordinate vector
8166858538eSMatthew G. Knepley 
817dce8aebaSBarry Smith   Level: intermediate
818dce8aebaSBarry Smith 
8196858538eSMatthew G. Knepley   Notes:
820dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
821dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
8226858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
8236858538eSMatthew G. Knepley 
8246858538eSMatthew G. Knepley   The vector c should be destroyed by the caller.
8256858538eSMatthew G. Knepley 
826dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
8276858538eSMatthew G. Knepley @*/
828d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
829d71ae5a4SJacob Faibussowitsch {
8306858538eSMatthew G. Knepley   PetscFunctionBegin;
8316858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8326858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
8336858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
8346858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
8356858538eSMatthew G. Knepley   dm->coordinates[1].xl = c;
8366858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
8373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8386858538eSMatthew G. Knepley }
8396858538eSMatthew G. Knepley 
840d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
841d71ae5a4SJacob Faibussowitsch {
8426858538eSMatthew G. Knepley   PetscFunctionBegin;
8436858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8446858538eSMatthew G. Knepley   PetscValidPointer(field, 2);
8456858538eSMatthew G. Knepley   if (!dm->coordinates[0].field) {
8466858538eSMatthew G. Knepley     if (dm->ops->createcoordinatefield) PetscCall((*dm->ops->createcoordinatefield)(dm, &dm->coordinates[0].field));
8476858538eSMatthew G. Knepley   }
8486858538eSMatthew G. Knepley   *field = dm->coordinates[0].field;
8493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8506858538eSMatthew G. Knepley }
8516858538eSMatthew G. Knepley 
852d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
853d71ae5a4SJacob Faibussowitsch {
8546858538eSMatthew G. Knepley   PetscFunctionBegin;
8556858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8566858538eSMatthew G. Knepley   if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2);
8576858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)field));
8586858538eSMatthew G. Knepley   PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
8596858538eSMatthew G. Knepley   dm->coordinates[0].field = field;
8603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8616858538eSMatthew G. Knepley }
8626858538eSMatthew G. Knepley 
8636858538eSMatthew G. Knepley /*@
864dce8aebaSBarry Smith   DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.
8656858538eSMatthew G. Knepley 
866*20f4b53cSBarry Smith   Not Collective
8676858538eSMatthew G. Knepley 
8686858538eSMatthew G. Knepley   Input Parameter:
869dce8aebaSBarry Smith . dm - the `DM`
8706858538eSMatthew G. Knepley 
8716858538eSMatthew G. Knepley   Output Parameters:
8726858538eSMatthew G. Knepley + lmin - local minimum coordinates (length coord dim, optional)
87335cb6cd3SPierre Jolivet - lmax - local maximum coordinates (length coord dim, optional)
8746858538eSMatthew G. Knepley 
8756858538eSMatthew G. Knepley   Level: beginner
8766858538eSMatthew G. Knepley 
877dce8aebaSBarry Smith   Note:
878dce8aebaSBarry Smith   If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.
8796858538eSMatthew G. Knepley 
880dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
8816858538eSMatthew G. Knepley @*/
882d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
883d71ae5a4SJacob Faibussowitsch {
8846858538eSMatthew G. Knepley   Vec       coords = NULL;
8856858538eSMatthew G. Knepley   PetscReal min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
8866858538eSMatthew G. Knepley   PetscReal max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
8876858538eSMatthew G. Knepley   PetscInt  cdim, i, j;
8886858538eSMatthew G. Knepley 
8896858538eSMatthew G. Knepley   PetscFunctionBegin;
8906858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8916858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
892bdf15c88SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
893bdf15c88SMatthew G. Knepley   if (coords) {
894bdf15c88SMatthew G. Knepley     const PetscScalar *local_coords;
895bdf15c88SMatthew G. Knepley     PetscInt           N, Ni;
896bdf15c88SMatthew G. Knepley 
897bdf15c88SMatthew G. Knepley     for (j = cdim; j < 3; ++j) {
898bdf15c88SMatthew G. Knepley       min[j] = 0;
899bdf15c88SMatthew G. Knepley       max[j] = 0;
900bdf15c88SMatthew G. Knepley     }
901bdf15c88SMatthew G. Knepley     PetscCall(VecGetArrayRead(coords, &local_coords));
902bdf15c88SMatthew G. Knepley     PetscCall(VecGetLocalSize(coords, &N));
903bdf15c88SMatthew G. Knepley     Ni = N / cdim;
904bdf15c88SMatthew G. Knepley     for (i = 0; i < Ni; ++i) {
905bdf15c88SMatthew G. Knepley       for (j = 0; j < cdim; ++j) {
906bdf15c88SMatthew G. Knepley         min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
907bdf15c88SMatthew G. Knepley         max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
908bdf15c88SMatthew G. Knepley       }
909bdf15c88SMatthew G. Knepley     }
910bdf15c88SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(coords, &local_coords));
911bdf15c88SMatthew G. Knepley     PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
9126858538eSMatthew G. Knepley     if (coords) {
9136858538eSMatthew G. Knepley       PetscCall(VecGetArrayRead(coords, &local_coords));
9146858538eSMatthew G. Knepley       PetscCall(VecGetLocalSize(coords, &N));
9156858538eSMatthew G. Knepley       Ni = N / cdim;
9166858538eSMatthew G. Knepley       for (i = 0; i < Ni; ++i) {
917bdf15c88SMatthew G. Knepley         for (j = 0; j < cdim; ++j) {
918bdf15c88SMatthew G. Knepley           min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
919bdf15c88SMatthew G. Knepley           max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
9206858538eSMatthew G. Knepley         }
9216858538eSMatthew G. Knepley       }
9226858538eSMatthew G. Knepley       PetscCall(VecRestoreArrayRead(coords, &local_coords));
923bdf15c88SMatthew G. Knepley     }
9246858538eSMatthew G. Knepley   } else {
9256858538eSMatthew G. Knepley     PetscBool isda;
9266858538eSMatthew G. Knepley 
9276858538eSMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &isda));
9286858538eSMatthew G. Knepley     if (isda) PetscCall(DMGetLocalBoundingIndices_DMDA(dm, min, max));
9296858538eSMatthew G. Knepley   }
9306858538eSMatthew G. Knepley   if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
9316858538eSMatthew G. Knepley   if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
9323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9336858538eSMatthew G. Knepley }
9346858538eSMatthew G. Knepley 
9356858538eSMatthew G. Knepley /*@
936dce8aebaSBarry Smith   DMGetBoundingBox - Returns the global bounding box for the `DM`.
9376858538eSMatthew G. Knepley 
9386858538eSMatthew G. Knepley   Collective
9396858538eSMatthew G. Knepley 
9406858538eSMatthew G. Knepley   Input Parameter:
941dce8aebaSBarry Smith . dm - the `DM`
9426858538eSMatthew G. Knepley 
9436858538eSMatthew G. Knepley   Output Parameters:
9446858538eSMatthew G. Knepley + gmin - global minimum coordinates (length coord dim, optional)
94535cb6cd3SPierre Jolivet - gmax - global maximum coordinates (length coord dim, optional)
9466858538eSMatthew G. Knepley 
9476858538eSMatthew G. Knepley   Level: beginner
9486858538eSMatthew G. Knepley 
949dce8aebaSBarry Smith .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
9506858538eSMatthew G. Knepley @*/
951d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
952d71ae5a4SJacob Faibussowitsch {
9536858538eSMatthew G. Knepley   PetscReal   lmin[3], lmax[3];
9546858538eSMatthew G. Knepley   PetscInt    cdim;
9556858538eSMatthew G. Knepley   PetscMPIInt count;
9566858538eSMatthew G. Knepley 
9576858538eSMatthew G. Knepley   PetscFunctionBegin;
9586858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9596858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
9606858538eSMatthew G. Knepley   PetscCall(PetscMPIIntCast(cdim, &count));
9616858538eSMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
9626858538eSMatthew G. Knepley   if (gmin) PetscCall(MPIU_Allreduce(lmin, gmin, count, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
9636858538eSMatthew G. Knepley   if (gmax) PetscCall(MPIU_Allreduce(lmax, gmax, count, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
9643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9656858538eSMatthew G. Knepley }
9666858538eSMatthew G. Knepley 
96790612307SJed Brown 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[])
96890612307SJed Brown {
96990612307SJed Brown   for (PetscInt i = 0; i < dim; i++) xnew[i] = x[i];
97090612307SJed Brown }
97190612307SJed Brown 
9726858538eSMatthew G. Knepley /*@
9736858538eSMatthew G. Knepley   DMProjectCoordinates - Project coordinates to a different space
9746858538eSMatthew G. Knepley 
9756858538eSMatthew G. Knepley   Input Parameters:
976dce8aebaSBarry Smith + dm      - The `DM` object
9776858538eSMatthew G. Knepley - disc    - The new coordinate discretization or NULL to ensure a coordinate discretization exists
9786858538eSMatthew G. Knepley 
9796858538eSMatthew G. Knepley   Level: intermediate
9806858538eSMatthew G. Knepley 
981dce8aebaSBarry Smith   Notes:
982dce8aebaSBarry Smith   A `PetscFE` defines an approximation space using a `PetscSpace`, which represents the basis functions, and a `PetscDualSpace`, which defines the interpolation operation
983dce8aebaSBarry Smith   in the space.
984dce8aebaSBarry Smith 
985dce8aebaSBarry Smith   This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
986dce8aebaSBarry Smith   The coordinate projection is done on the continuous coordinates, and if possible, the discontinuous coordinates are also updated.
987dce8aebaSBarry Smith 
988dce8aebaSBarry Smith   Developer Note:
989dce8aebaSBarry Smith   With more effort, we could directly project the discontinuous coordinates also.
990dce8aebaSBarry Smith 
991dce8aebaSBarry Smith .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
9926858538eSMatthew G. Knepley @*/
993d71ae5a4SJacob Faibussowitsch PetscErrorCode DMProjectCoordinates(DM dm, PetscFE disc)
994d71ae5a4SJacob Faibussowitsch {
9956858538eSMatthew G. Knepley   PetscFE      discOld;
9966858538eSMatthew G. Knepley   PetscClassId classid;
9976858538eSMatthew G. Knepley   DM           cdmOld, cdmNew;
9986858538eSMatthew G. Knepley   Vec          coordsOld, coordsNew;
9996858538eSMatthew G. Knepley   PetscBool    same_space = PETSC_TRUE;
1000dd4c3f67SMatthew G. Knepley   const char  *prefix;
10016858538eSMatthew G. Knepley 
10026858538eSMatthew G. Knepley   PetscFunctionBegin;
10036858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10046858538eSMatthew G. Knepley   if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2);
10056858538eSMatthew G. Knepley 
10066858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdmOld));
10076858538eSMatthew G. Knepley   /* Check current discretization is compatible */
10086858538eSMatthew G. Knepley   PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
10096858538eSMatthew G. Knepley   PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
10106858538eSMatthew G. Knepley   if (classid != PETSCFE_CLASSID) {
10116858538eSMatthew G. Knepley     if (classid == PETSC_CONTAINER_CLASSID) {
10126858538eSMatthew G. Knepley       PetscFE        feLinear;
10136858538eSMatthew G. Knepley       DMPolytopeType ct;
1014dc431b0cSMatthew G. Knepley       PetscInt       dim, dE, cStart, cEnd, ctTmp;
10156858538eSMatthew G. Knepley 
10166858538eSMatthew G. Knepley       /* Assume linear vertex coordinates */
10176858538eSMatthew G. Knepley       PetscCall(DMGetDimension(dm, &dim));
10186858538eSMatthew G. Knepley       PetscCall(DMGetCoordinateDim(dm, &dE));
10196858538eSMatthew G. Knepley       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1020dc431b0cSMatthew G. Knepley       if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1021dc431b0cSMatthew G. Knepley       else ct = DM_POLYTOPE_UNKNOWN;
1022dc431b0cSMatthew G. Knepley       ctTmp = (PetscInt)ct;
1023dc431b0cSMatthew G. Knepley       PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &ctTmp, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1024dc431b0cSMatthew G. Knepley       ct = (DMPolytopeType)ctTmp;
1025dc431b0cSMatthew G. Knepley       PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, 1, -1, &feLinear));
10266858538eSMatthew G. Knepley       PetscCall(DMSetField(cdmOld, 0, NULL, (PetscObject)feLinear));
10276858538eSMatthew G. Knepley       PetscCall(PetscFEDestroy(&feLinear));
10286858538eSMatthew G. Knepley       PetscCall(DMCreateDS(cdmOld));
10296858538eSMatthew G. Knepley       PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
10306858538eSMatthew G. Knepley     } else {
10316858538eSMatthew G. Knepley       const char *discname;
10326858538eSMatthew G. Knepley 
10336858538eSMatthew G. Knepley       PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
10346858538eSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
10356858538eSMatthew G. Knepley     }
10366858538eSMatthew G. Knepley   }
10373ba16761SJacob Faibussowitsch   if (!disc) PetscFunctionReturn(PETSC_SUCCESS);
10386858538eSMatthew G. Knepley   { // Check if the new space is the same as the old modulo quadrature
10396858538eSMatthew G. Knepley     PetscDualSpace dsOld, ds;
10406858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
10416858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(disc, &ds));
10426858538eSMatthew G. Knepley     PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
10436858538eSMatthew G. Knepley   }
10446858538eSMatthew G. Knepley   /* Make a fresh clone of the coordinate DM */
10456858538eSMatthew G. Knepley   PetscCall(DMClone(cdmOld, &cdmNew));
1046dd4c3f67SMatthew G. Knepley   cdmNew->cloneOpts = PETSC_TRUE;
1047dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix));
1048dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix));
10496858538eSMatthew G. Knepley   PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
10506858538eSMatthew G. Knepley   PetscCall(DMCreateDS(cdmNew));
10516725e60dSJed Brown   if (cdmOld->periodic.setup) {
10526725e60dSJed Brown     cdmNew->periodic.setup = cdmOld->periodic.setup;
10536725e60dSJed Brown     PetscCall(cdmNew->periodic.setup(cdmNew));
10546725e60dSJed Brown   }
1055dd4c3f67SMatthew G. Knepley   if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew));
10566858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coordsOld));
10576858538eSMatthew G. Knepley   PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
1058f0ac6e35SMatthew G. Knepley   if (same_space) {
1059f0ac6e35SMatthew G. Knepley     // Need to copy so that the new vector has the right dm
1060f0ac6e35SMatthew G. Knepley     PetscCall(VecCopy(coordsOld, coordsNew));
106190612307SJed Brown   } else { // Project the coordinate vector from old to new space
106290612307SJed Brown     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};
106390612307SJed Brown     // We can't call DMProjectField directly because it depends on KSP for DMGlobalToLocalSolve(), but we can use the core strategy
106490612307SJed Brown     Vec X_new_loc;
106590612307SJed Brown     PetscCall(DMCreateLocalVector(cdmNew, &X_new_loc));
106690612307SJed Brown     PetscCall(DMSetCoordinateDM(cdmNew, cdmOld));
106790612307SJed Brown     // See DMPlexRemapGeometry() for a similar pattern handling the coordinate field
106890612307SJed Brown     DMField cf;
106990612307SJed Brown     PetscCall(DMGetCoordinateField(dm, &cf));
107090612307SJed Brown     cdmNew->coordinates[0].field = cf;
107190612307SJed Brown     PetscCall(DMProjectFieldLocal(cdmNew, 0.0, NULL, funcs, INSERT_VALUES, X_new_loc));
107290612307SJed Brown     cdmNew->coordinates[0].field = NULL;
107390612307SJed Brown     PetscCall(DMSetCoordinateDM(cdmNew, NULL));
107490612307SJed Brown     PetscCall(DMLocalToGlobal(cdmNew, X_new_loc, INSERT_VALUES, coordsNew));
107590612307SJed Brown     PetscCall(VecDestroy(&X_new_loc));
10766858538eSMatthew G. Knepley   }
10776858538eSMatthew G. Knepley   /* Set new coordinate structures */
10786858538eSMatthew G. Knepley   PetscCall(DMSetCoordinateField(dm, NULL));
10796858538eSMatthew G. Knepley   PetscCall(DMSetCoordinateDM(dm, cdmNew));
10806858538eSMatthew G. Knepley   PetscCall(DMSetCoordinates(dm, coordsNew));
10816858538eSMatthew G. Knepley   PetscCall(VecDestroy(&coordsNew));
10826858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cdmNew));
10833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10846858538eSMatthew G. Knepley }
10856858538eSMatthew G. Knepley 
10866858538eSMatthew G. Knepley /*@
1087dce8aebaSBarry Smith   DMLocatePoints - Locate the points in v in the mesh and return a `PetscSF` of the containing cells
10886858538eSMatthew G. Knepley 
1089*20f4b53cSBarry Smith   Collective
10906858538eSMatthew G. Knepley 
10916858538eSMatthew G. Knepley   Input Parameters:
1092dce8aebaSBarry Smith + dm - The `DM`
1093dce8aebaSBarry Smith - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`
10946858538eSMatthew G. Knepley 
10956858538eSMatthew G. Knepley   Input/Output Parameters:
1096dce8aebaSBarry Smith + v - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
1097*20f4b53cSBarry Smith - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point;
1098dce8aebaSBarry Smith            on output, the `PetscSF` containing the ranks and local indices of the containing points
10996858538eSMatthew G. Knepley 
11006858538eSMatthew G. Knepley   Level: developer
11016858538eSMatthew G. Knepley 
11026858538eSMatthew G. Knepley   Notes:
11034a8fad2eSMatthew G. Knepley   To do a search of the local cells of the mesh, v should have `PETSC_COMM_SELF` as its communicator.
1104*20f4b53cSBarry Smith   To do a search of all the cells in the distributed mesh, `v` should have the same communicator as `dm`.
11056858538eSMatthew G. Knepley 
1106d8206211SMatthew G. Knepley   Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.
1107d8206211SMatthew G. Knepley 
1108*20f4b53cSBarry Smith   If *cellSF is `NULL` on input, a `PetscSF` will be created.
1109*20f4b53cSBarry Smith   If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.
11106858538eSMatthew G. Knepley 
11116858538eSMatthew G. Knepley   An array that maps each point to its containing cell can be obtained with
1112dce8aebaSBarry Smith .vb
1113dce8aebaSBarry Smith     const PetscSFNode *cells;
1114dce8aebaSBarry Smith     PetscInt           nFound;
1115dce8aebaSBarry Smith     const PetscInt    *found;
11166858538eSMatthew G. Knepley 
1117dce8aebaSBarry Smith     PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1118dce8aebaSBarry Smith .ve
11196858538eSMatthew G. Knepley 
11204a8fad2eSMatthew G. Knepley   Where cells[i].rank is the rank of the process owning the cell containing point found[i] (or i if found == NULL), and cells[i].index is
1121*20f4b53cSBarry Smith   the index of the cell in its rank's 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.
11226858538eSMatthew G. Knepley 
1123dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
11246858538eSMatthew G. Knepley @*/
1125d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1126d71ae5a4SJacob Faibussowitsch {
11276858538eSMatthew G. Knepley   PetscFunctionBegin;
11286858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11296858538eSMatthew G. Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 2);
11306858538eSMatthew G. Knepley   PetscValidPointer(cellSF, 4);
11316858538eSMatthew G. Knepley   if (*cellSF) {
11326858538eSMatthew G. Knepley     PetscMPIInt result;
11336858538eSMatthew G. Knepley 
11346858538eSMatthew G. Knepley     PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4);
11356858538eSMatthew G. Knepley     PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
11366858538eSMatthew 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");
11376858538eSMatthew G. Knepley   } else {
11386858538eSMatthew G. Knepley     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
11396858538eSMatthew G. Knepley   }
11406858538eSMatthew G. Knepley   PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1141dbbe0bcdSBarry Smith   PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
11426858538eSMatthew G. Knepley   PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
11433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11446858538eSMatthew G. Knepley }
1145