xref: /petsc/src/dm/interface/dmcoordinates.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
16858538eSMatthew G. Knepley #include <petsc/private/dmimpl.h> /*I      "petscdm.h"          I*/
26858538eSMatthew G. Knepley 
3e44f6aebSMatthew G. Knepley #include <petscdmplex.h> /* For DMCreateAffineCoordinates_Internal() */
46858538eSMatthew G. Knepley #include <petscsf.h>     /* For DMLocatePoints() */
56858538eSMatthew G. Knepley 
DMRestrictHook_Coordinates(DM dm,DM dmc,PetscCtx ctx)6*2a8381b2SBarry Smith PetscErrorCode DMRestrictHook_Coordinates(DM dm, DM dmc, PetscCtx ctx)
7d71ae5a4SJacob Faibussowitsch {
86858538eSMatthew G. Knepley   DM  dm_coord, dmc_coord;
96858538eSMatthew G. Knepley   Vec coords, ccoords;
106858538eSMatthew G. Knepley   Mat inject;
114d86920dSPierre Jolivet 
126858538eSMatthew G. Knepley   PetscFunctionBegin;
136858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
146858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dmc, &dmc_coord));
156858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
166858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dmc, &ccoords));
176858538eSMatthew G. Knepley   if (coords && !ccoords) {
186858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(dmc_coord, &ccoords));
196858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
206858538eSMatthew G. Knepley     PetscCall(DMCreateInjection(dmc_coord, dm_coord, &inject));
216858538eSMatthew G. Knepley     PetscCall(MatRestrict(inject, coords, ccoords));
226858538eSMatthew G. Knepley     PetscCall(MatDestroy(&inject));
236858538eSMatthew G. Knepley     PetscCall(DMSetCoordinates(dmc, ccoords));
246858538eSMatthew G. Knepley     PetscCall(VecDestroy(&ccoords));
256858538eSMatthew G. Knepley   }
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
276858538eSMatthew G. Knepley }
286858538eSMatthew G. Knepley 
DMSubDomainHook_Coordinates(DM dm,DM subdm,PetscCtx ctx)29*2a8381b2SBarry Smith static PetscErrorCode DMSubDomainHook_Coordinates(DM dm, DM subdm, PetscCtx ctx)
30d71ae5a4SJacob Faibussowitsch {
316858538eSMatthew G. Knepley   DM          dm_coord, subdm_coord;
326858538eSMatthew G. Knepley   Vec         coords, ccoords, clcoords;
336858538eSMatthew G. Knepley   VecScatter *scat_i, *scat_g;
344d86920dSPierre Jolivet 
356858538eSMatthew G. Knepley   PetscFunctionBegin;
366858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
376858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(subdm, &subdm_coord));
386858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(dm, &coords));
396858538eSMatthew G. Knepley   PetscCall(DMGetCoordinates(subdm, &ccoords));
406858538eSMatthew G. Knepley   if (coords && !ccoords) {
416858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(subdm_coord, &ccoords));
426858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)ccoords, "coordinates"));
436858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(subdm_coord, &clcoords));
446858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)clcoords, "coordinates"));
456858538eSMatthew G. Knepley     PetscCall(DMCreateDomainDecompositionScatters(dm_coord, 1, &subdm_coord, NULL, &scat_i, &scat_g));
466858538eSMatthew G. Knepley     PetscCall(VecScatterBegin(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
476858538eSMatthew G. Knepley     PetscCall(VecScatterEnd(scat_i[0], coords, ccoords, INSERT_VALUES, SCATTER_FORWARD));
486858538eSMatthew G. Knepley     PetscCall(VecScatterBegin(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
496858538eSMatthew G. Knepley     PetscCall(VecScatterEnd(scat_g[0], coords, clcoords, INSERT_VALUES, SCATTER_FORWARD));
506858538eSMatthew G. Knepley     PetscCall(DMSetCoordinates(subdm, ccoords));
516858538eSMatthew G. Knepley     PetscCall(DMSetCoordinatesLocal(subdm, clcoords));
526858538eSMatthew G. Knepley     PetscCall(VecScatterDestroy(&scat_i[0]));
536858538eSMatthew G. Knepley     PetscCall(VecScatterDestroy(&scat_g[0]));
546858538eSMatthew G. Knepley     PetscCall(VecDestroy(&ccoords));
556858538eSMatthew G. Knepley     PetscCall(VecDestroy(&clcoords));
566858538eSMatthew G. Knepley     PetscCall(PetscFree(scat_i));
576858538eSMatthew G. Knepley     PetscCall(PetscFree(scat_g));
586858538eSMatthew G. Knepley   }
593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
606858538eSMatthew G. Knepley }
616858538eSMatthew G. Knepley 
626858538eSMatthew G. Knepley /*@
63dce8aebaSBarry Smith   DMGetCoordinateDM - Gets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
646858538eSMatthew G. Knepley 
6520f4b53cSBarry Smith   Collective
666858538eSMatthew G. Knepley 
676858538eSMatthew G. Knepley   Input Parameter:
68dce8aebaSBarry Smith . dm - the `DM`
696858538eSMatthew G. Knepley 
706858538eSMatthew G. Knepley   Output Parameter:
71dce8aebaSBarry Smith . cdm - coordinate `DM`
726858538eSMatthew G. Knepley 
736858538eSMatthew G. Knepley   Level: intermediate
746858538eSMatthew G. Knepley 
75dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDM()`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGSetCellCoordinateDM()`,
7660225df5SJacob Faibussowitsch 
776858538eSMatthew G. Knepley @*/
DMGetCoordinateDM(DM dm,DM * cdm)78d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDM(DM dm, DM *cdm)
79d71ae5a4SJacob Faibussowitsch {
806858538eSMatthew G. Knepley   PetscFunctionBegin;
816858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
824f572ea9SToby Isaac   PetscAssertPointer(cdm, 2);
836858538eSMatthew G. Knepley   if (!dm->coordinates[0].dm) {
846858538eSMatthew G. Knepley     DM cdm;
856858538eSMatthew G. Knepley 
86dbbe0bcdSBarry Smith     PetscUseTypeMethod(dm, createcoordinatedm, &cdm);
876858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)cdm, "coordinateDM"));
886858538eSMatthew G. Knepley     /* Just in case the DM sets the coordinate DM when creating it (DMP4est can do this, because it may not setup
896858538eSMatthew G. Knepley      * until the call to CreateCoordinateDM) */
906858538eSMatthew G. Knepley     PetscCall(DMDestroy(&dm->coordinates[0].dm));
916858538eSMatthew G. Knepley     dm->coordinates[0].dm = cdm;
926858538eSMatthew G. Knepley   }
936858538eSMatthew G. Knepley   *cdm = dm->coordinates[0].dm;
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
956858538eSMatthew G. Knepley }
966858538eSMatthew G. Knepley 
976858538eSMatthew G. Knepley /*@
98dce8aebaSBarry Smith   DMSetCoordinateDM - Sets the `DM` that prescribes coordinate layout and scatters between global and local coordinates
996858538eSMatthew G. Knepley 
10020f4b53cSBarry Smith   Logically Collective
1016858538eSMatthew G. Knepley 
1026858538eSMatthew G. Knepley   Input Parameters:
103dce8aebaSBarry Smith + dm  - the `DM`
104dce8aebaSBarry Smith - cdm - coordinate `DM`
1056858538eSMatthew G. Knepley 
1066858538eSMatthew G. Knepley   Level: intermediate
1076858538eSMatthew G. Knepley 
108dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMSetCoordinates()`, `DMGetCellCoordinateDM()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`,
109dce8aebaSBarry Smith           `DMGSetCellCoordinateDM()`
1106858538eSMatthew G. Knepley @*/
DMSetCoordinateDM(DM dm,DM cdm)111d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDM(DM dm, DM cdm)
112d71ae5a4SJacob Faibussowitsch {
1136858538eSMatthew G. Knepley   PetscFunctionBegin;
1146858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11590612307SJed Brown   if (cdm) PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
1166858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)cdm));
1176858538eSMatthew G. Knepley   PetscCall(DMDestroy(&dm->coordinates[0].dm));
1186858538eSMatthew G. Knepley   dm->coordinates[0].dm = cdm;
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1206858538eSMatthew G. Knepley }
1216858538eSMatthew G. Knepley 
1226858538eSMatthew G. Knepley /*@
123dce8aebaSBarry Smith   DMGetCellCoordinateDM - Gets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1246858538eSMatthew G. Knepley 
12520f4b53cSBarry Smith   Collective
1266858538eSMatthew G. Knepley 
1276858538eSMatthew G. Knepley   Input Parameter:
128dce8aebaSBarry Smith . dm - the `DM`
1296858538eSMatthew G. Knepley 
1306858538eSMatthew G. Knepley   Output Parameter:
1310b3275a6SBarry Smith . cdm - cellwise coordinate `DM`, or `NULL` if they are not defined
1326858538eSMatthew G. Knepley 
1336858538eSMatthew G. Knepley   Level: intermediate
1346858538eSMatthew G. Knepley 
135dce8aebaSBarry Smith   Note:
136dce8aebaSBarry Smith   Call `DMLocalizeCoordinates()` to automatically create cellwise coordinates for periodic geometries.
137dce8aebaSBarry Smith 
138dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
139dce8aebaSBarry Smith           `DMLocalizeCoordinates()`, `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
1406858538eSMatthew G. Knepley @*/
DMGetCellCoordinateDM(DM dm,DM * cdm)141d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateDM(DM dm, DM *cdm)
142d71ae5a4SJacob Faibussowitsch {
1436858538eSMatthew G. Knepley   PetscFunctionBegin;
1446858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1454f572ea9SToby Isaac   PetscAssertPointer(cdm, 2);
1466858538eSMatthew G. Knepley   *cdm = dm->coordinates[1].dm;
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1486858538eSMatthew G. Knepley }
1496858538eSMatthew G. Knepley 
1506858538eSMatthew G. Knepley /*@
151dce8aebaSBarry Smith   DMSetCellCoordinateDM - Sets the `DM` that prescribes cellwise coordinate layout and scatters between global and local cellwise coordinates
1526858538eSMatthew G. Knepley 
15320f4b53cSBarry Smith   Logically Collective
1546858538eSMatthew G. Knepley 
1556858538eSMatthew G. Knepley   Input Parameters:
156dce8aebaSBarry Smith + dm  - the `DM`
157dce8aebaSBarry Smith - cdm - cellwise coordinate `DM`
1586858538eSMatthew G. Knepley 
1596858538eSMatthew G. Knepley   Level: intermediate
1606858538eSMatthew G. Knepley 
161dce8aebaSBarry Smith   Note:
16235cb6cd3SPierre Jolivet   As opposed to `DMSetCoordinateDM()` these coordinates are useful for discontinuous Galerkin methods since they support coordinate fields that are discontinuous at cell boundaries.
163dce8aebaSBarry Smith 
164dce8aebaSBarry Smith .seealso: `DMGetCellCoordinateDM()`, `DMSetCellCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`,
165dce8aebaSBarry Smith           `DMSetCoordinateDM()`, `DMGetCoordinateDM()`
1666858538eSMatthew G. Knepley @*/
DMSetCellCoordinateDM(DM dm,DM cdm)167d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateDM(DM dm, DM cdm)
168d71ae5a4SJacob Faibussowitsch {
1696858538eSMatthew G. Knepley   PetscInt dim;
1706858538eSMatthew G. Knepley 
1716858538eSMatthew G. Knepley   PetscFunctionBegin;
1726858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1736858538eSMatthew G. Knepley   if (cdm) {
1746858538eSMatthew G. Knepley     PetscValidHeaderSpecific(cdm, DM_CLASSID, 2);
1756858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDim(dm, &dim));
1766858538eSMatthew G. Knepley     dm->coordinates[1].dim = dim;
1776858538eSMatthew G. Knepley   }
1786858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)cdm));
1796858538eSMatthew G. Knepley   PetscCall(DMDestroy(&dm->coordinates[1].dm));
1806858538eSMatthew G. Knepley   dm->coordinates[1].dm = cdm;
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1826858538eSMatthew G. Knepley }
1836858538eSMatthew G. Knepley 
1846858538eSMatthew G. Knepley /*@
1850b3275a6SBarry Smith   DMGetCoordinateDim - Retrieve the dimension of the embedding space for coordinate values. For example a mesh on the surface of a sphere would have a 3 dimensional embedding space
1866858538eSMatthew G. Knepley 
1876858538eSMatthew G. Knepley   Not Collective
1886858538eSMatthew G. Knepley 
1896858538eSMatthew G. Knepley   Input Parameter:
190dce8aebaSBarry Smith . dm - The `DM` object
1916858538eSMatthew G. Knepley 
1926858538eSMatthew G. Knepley   Output Parameter:
1936858538eSMatthew G. Knepley . dim - The embedding dimension
1946858538eSMatthew G. Knepley 
1956858538eSMatthew G. Knepley   Level: intermediate
1966858538eSMatthew G. Knepley 
197dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
1986858538eSMatthew G. Knepley @*/
DMGetCoordinateDim(DM dm,PetscInt * dim)199d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateDim(DM dm, PetscInt *dim)
200d71ae5a4SJacob Faibussowitsch {
2016858538eSMatthew G. Knepley   PetscFunctionBegin;
2026858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2034f572ea9SToby Isaac   PetscAssertPointer(dim, 2);
2046858538eSMatthew G. Knepley   if (dm->coordinates[0].dim == PETSC_DEFAULT) dm->coordinates[0].dim = dm->dim;
2056858538eSMatthew G. Knepley   *dim = dm->coordinates[0].dim;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2076858538eSMatthew G. Knepley }
2086858538eSMatthew G. Knepley 
2096858538eSMatthew G. Knepley /*@
2106858538eSMatthew G. Knepley   DMSetCoordinateDim - Set the dimension of the embedding space for coordinate values.
2116858538eSMatthew G. Knepley 
2126858538eSMatthew G. Knepley   Not Collective
2136858538eSMatthew G. Knepley 
2146858538eSMatthew G. Knepley   Input Parameters:
215dce8aebaSBarry Smith + dm  - The `DM` object
2166858538eSMatthew G. Knepley - dim - The embedding dimension
2176858538eSMatthew G. Knepley 
2186858538eSMatthew G. Knepley   Level: intermediate
2196858538eSMatthew G. Knepley 
220dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2216858538eSMatthew G. Knepley @*/
DMSetCoordinateDim(DM dm,PetscInt dim)222d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateDim(DM dm, PetscInt dim)
223d71ae5a4SJacob Faibussowitsch {
2246858538eSMatthew G. Knepley   PetscDS  ds;
2256858538eSMatthew G. Knepley   PetscInt Nds, n;
2266858538eSMatthew G. Knepley 
2276858538eSMatthew G. Knepley   PetscFunctionBegin;
2286858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2296858538eSMatthew G. Knepley   dm->coordinates[0].dim = dim;
2306858538eSMatthew G. Knepley   if (dm->dim >= 0) {
2316858538eSMatthew G. Knepley     PetscCall(DMGetNumDS(dm, &Nds));
2326858538eSMatthew G. Knepley     for (n = 0; n < Nds; ++n) {
23307218a29SMatthew G. Knepley       PetscCall(DMGetRegionNumDS(dm, n, NULL, NULL, &ds, NULL));
2346858538eSMatthew G. Knepley       PetscCall(PetscDSSetCoordinateDimension(ds, dim));
2356858538eSMatthew G. Knepley     }
2366858538eSMatthew G. Knepley   }
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2386858538eSMatthew G. Knepley }
2396858538eSMatthew G. Knepley 
2406858538eSMatthew G. Knepley /*@
2410b3275a6SBarry Smith   DMGetCoordinateSection - Retrieve the `PetscSection` of coordinate values over the mesh.
2426858538eSMatthew G. Knepley 
24320f4b53cSBarry Smith   Collective
2446858538eSMatthew G. Knepley 
2456858538eSMatthew G. Knepley   Input Parameter:
246dce8aebaSBarry Smith . dm - The `DM` object
2476858538eSMatthew G. Knepley 
2486858538eSMatthew G. Knepley   Output Parameter:
249dce8aebaSBarry Smith . section - The `PetscSection` object
2506858538eSMatthew G. Knepley 
2516858538eSMatthew G. Knepley   Level: intermediate
2526858538eSMatthew G. Knepley 
253dce8aebaSBarry Smith   Note:
254dce8aebaSBarry Smith   This just retrieves the local section from the coordinate `DM`. In other words,
255dce8aebaSBarry Smith .vb
256dce8aebaSBarry Smith   DMGetCoordinateDM(dm, &cdm);
257dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
258dce8aebaSBarry Smith .ve
259dce8aebaSBarry Smith 
2600b3275a6SBarry Smith .seealso: `DM`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2616858538eSMatthew G. Knepley @*/
DMGetCoordinateSection(DM dm,PetscSection * section)262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateSection(DM dm, PetscSection *section)
263d71ae5a4SJacob Faibussowitsch {
2646858538eSMatthew G. Knepley   DM cdm;
2656858538eSMatthew G. Knepley 
2666858538eSMatthew G. Knepley   PetscFunctionBegin;
2676858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2684f572ea9SToby Isaac   PetscAssertPointer(section, 2);
2696858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2706858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, section));
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2726858538eSMatthew G. Knepley }
2736858538eSMatthew G. Knepley 
2746858538eSMatthew G. Knepley /*@
2750b3275a6SBarry Smith   DMSetCoordinateSection - Set the `PetscSection` of coordinate values over the mesh.
2766858538eSMatthew G. Knepley 
2776858538eSMatthew G. Knepley   Not Collective
2786858538eSMatthew G. Knepley 
2796858538eSMatthew G. Knepley   Input Parameters:
280dce8aebaSBarry Smith + dm      - The `DM` object
281dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
282dce8aebaSBarry Smith - section - The `PetscSection` object
2836858538eSMatthew G. Knepley 
2846858538eSMatthew G. Knepley   Level: intermediate
2856858538eSMatthew G. Knepley 
286dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMGetCoordinateSection()`, `DMGetLocalSection()`, `DMSetLocalSection()`
2876858538eSMatthew G. Knepley @*/
DMSetCoordinateSection(DM dm,PetscInt dim,PetscSection section)288d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateSection(DM dm, PetscInt dim, PetscSection section)
289d71ae5a4SJacob Faibussowitsch {
2906858538eSMatthew G. Knepley   DM cdm;
2916858538eSMatthew G. Knepley 
2926858538eSMatthew G. Knepley   PetscFunctionBegin;
2936858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2946858538eSMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
2956858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
2966858538eSMatthew G. Knepley   PetscCall(DMSetLocalSection(cdm, section));
2976858538eSMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
2986858538eSMatthew G. Knepley     PetscInt d = PETSC_DEFAULT;
2996858538eSMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
3006858538eSMatthew G. Knepley 
3016858538eSMatthew G. Knepley     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3026858538eSMatthew G. Knepley     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
3036858538eSMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
3046858538eSMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
3056858538eSMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
3066858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(section, v, &dd));
3079371c9d4SSatish Balay       if (dd) {
3089371c9d4SSatish Balay         d = dd;
3099371c9d4SSatish Balay         break;
3109371c9d4SSatish Balay       }
3116858538eSMatthew G. Knepley     }
3126858538eSMatthew G. Knepley     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
3134c712d99Sksagiyam   } else {
3144c712d99Sksagiyam     PetscCall(DMSetCoordinateDim(dm, dim));
3156858538eSMatthew G. Knepley   }
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3176858538eSMatthew G. Knepley }
3186858538eSMatthew G. Knepley 
3196858538eSMatthew G. Knepley /*@
3200b3275a6SBarry Smith   DMGetCellCoordinateSection - Retrieve the `PetscSection` of cellwise coordinate values over the mesh.
3216858538eSMatthew G. Knepley 
32220f4b53cSBarry Smith   Collective
3236858538eSMatthew G. Knepley 
3246858538eSMatthew G. Knepley   Input Parameter:
325dce8aebaSBarry Smith . dm - The `DM` object
3266858538eSMatthew G. Knepley 
3276858538eSMatthew G. Knepley   Output Parameter:
32820f4b53cSBarry Smith . section - The `PetscSection` object, or `NULL` if no cellwise coordinates are defined
3296858538eSMatthew G. Knepley 
3306858538eSMatthew G. Knepley   Level: intermediate
3316858538eSMatthew G. Knepley 
332dce8aebaSBarry Smith   Note:
333dce8aebaSBarry Smith   This just retrieves the local section from the cell coordinate `DM`. In other words,
334dce8aebaSBarry Smith .vb
335dce8aebaSBarry Smith   DMGetCellCoordinateDM(dm, &cdm);
336dce8aebaSBarry Smith   DMGetLocalSection(cdm, &section);
337dce8aebaSBarry Smith .ve
338dce8aebaSBarry Smith 
339dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateSection()`, `DMSetCellCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
3406858538eSMatthew G. Knepley @*/
DMGetCellCoordinateSection(DM dm,PetscSection * section)341d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinateSection(DM dm, PetscSection *section)
342d71ae5a4SJacob Faibussowitsch {
3436858538eSMatthew G. Knepley   DM cdm;
3446858538eSMatthew G. Knepley 
3456858538eSMatthew G. Knepley   PetscFunctionBegin;
3466858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3474f572ea9SToby Isaac   PetscAssertPointer(section, 2);
3486858538eSMatthew G. Knepley   *section = NULL;
3496858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3506858538eSMatthew G. Knepley   if (cdm) PetscCall(DMGetLocalSection(cdm, section));
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3526858538eSMatthew G. Knepley }
3536858538eSMatthew G. Knepley 
3546858538eSMatthew G. Knepley /*@
3550b3275a6SBarry Smith   DMSetCellCoordinateSection - Set the `PetscSection` of cellwise coordinate values over the mesh.
3566858538eSMatthew G. Knepley 
3576858538eSMatthew G. Knepley   Not Collective
3586858538eSMatthew G. Knepley 
3596858538eSMatthew G. Knepley   Input Parameters:
360dce8aebaSBarry Smith + dm      - The `DM` object
361dce8aebaSBarry Smith . dim     - The embedding dimension, or `PETSC_DETERMINE`
362dce8aebaSBarry Smith - section - The `PetscSection` object for a cellwise layout
3636858538eSMatthew G. Knepley 
3646858538eSMatthew G. Knepley   Level: intermediate
3656858538eSMatthew G. Knepley 
366dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinateDim()`, `DMSetCoordinateSection()`, `DMGetCellCoordinateSection()`, `DMGetCoordinateSection()`, `DMGetCellCoordinateDM()`, `DMGetLocalSection()`, `DMSetLocalSection()`
3676858538eSMatthew G. Knepley @*/
DMSetCellCoordinateSection(DM dm,PetscInt dim,PetscSection section)368d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinateSection(DM dm, PetscInt dim, PetscSection section)
369d71ae5a4SJacob Faibussowitsch {
3706858538eSMatthew G. Knepley   DM cdm;
3716858538eSMatthew G. Knepley 
3726858538eSMatthew G. Knepley   PetscFunctionBegin;
3736858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3746858538eSMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
3756858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3766858538eSMatthew G. Knepley   PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No DM defined for cellwise coordinates");
3776858538eSMatthew G. Knepley   PetscCall(DMSetLocalSection(cdm, section));
3786858538eSMatthew G. Knepley   if (dim == PETSC_DETERMINE) {
3796858538eSMatthew G. Knepley     PetscInt d = PETSC_DEFAULT;
3806858538eSMatthew G. Knepley     PetscInt pStart, pEnd, vStart, vEnd, v, dd;
3816858538eSMatthew G. Knepley 
3826858538eSMatthew G. Knepley     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
3836858538eSMatthew G. Knepley     PetscCall(DMGetDimPoints(dm, 0, &vStart, &vEnd));
3846858538eSMatthew G. Knepley     pStart = PetscMax(vStart, pStart);
3856858538eSMatthew G. Knepley     pEnd   = PetscMin(vEnd, pEnd);
3866858538eSMatthew G. Knepley     for (v = pStart; v < pEnd; ++v) {
3876858538eSMatthew G. Knepley       PetscCall(PetscSectionGetDof(section, v, &dd));
3889371c9d4SSatish Balay       if (dd) {
3899371c9d4SSatish Balay         d = dd;
3909371c9d4SSatish Balay         break;
3919371c9d4SSatish Balay       }
3926858538eSMatthew G. Knepley     }
3936858538eSMatthew G. Knepley     if (d >= 0) PetscCall(DMSetCoordinateDim(dm, d));
3944c712d99Sksagiyam   } else {
3954c712d99Sksagiyam     PetscCall(DMSetCoordinateDim(dm, dim));
3966858538eSMatthew G. Knepley   }
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3986858538eSMatthew G. Knepley }
3996858538eSMatthew G. Knepley 
4006858538eSMatthew G. Knepley /*@
401dce8aebaSBarry Smith   DMGetCoordinates - Gets a global vector with the coordinates associated with the `DM`.
4026858538eSMatthew G. Knepley 
403356ea6bcSBarry Smith   Collective if the global vector with coordinates has not been set yet but the local vector with coordinates has been set
4046858538eSMatthew G. Knepley 
4056858538eSMatthew G. Knepley   Input Parameter:
406dce8aebaSBarry Smith . dm - the `DM`
4076858538eSMatthew G. Knepley 
4086858538eSMatthew G. Knepley   Output Parameter:
4096858538eSMatthew G. Knepley . c - global coordinate vector
4106858538eSMatthew G. Knepley 
411dce8aebaSBarry Smith   Level: intermediate
412dce8aebaSBarry Smith 
413dce8aebaSBarry Smith   Notes:
414dce8aebaSBarry Smith   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
4150b3275a6SBarry Smith   destroyed `c` will no longer be valid.
4166858538eSMatthew G. Knepley 
417356ea6bcSBarry Smith   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates), see `DMGetCoordinatesLocal()`.
4186858538eSMatthew G. Knepley 
419dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
4206858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
4216858538eSMatthew G. Knepley 
422356ea6bcSBarry Smith   Does not work for `DMSTAG`
423356ea6bcSBarry Smith 
424dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
4256858538eSMatthew G. Knepley @*/
DMGetCoordinates(DM dm,Vec * c)426d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinates(DM dm, Vec *c)
427d71ae5a4SJacob Faibussowitsch {
4286858538eSMatthew G. Knepley   PetscFunctionBegin;
4296858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4304f572ea9SToby Isaac   PetscAssertPointer(c, 2);
4316858538eSMatthew G. Knepley   if (!dm->coordinates[0].x && dm->coordinates[0].xl) {
4326858538eSMatthew G. Knepley     DM cdm = NULL;
4336858538eSMatthew G. Knepley 
4346858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
4356858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[0].x));
4366858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].x, "coordinates"));
4376858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
4386858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[0].xl, INSERT_VALUES, dm->coordinates[0].x));
4396858538eSMatthew G. Knepley   }
4406858538eSMatthew G. Knepley   *c = dm->coordinates[0].x;
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4426858538eSMatthew G. Knepley }
4436858538eSMatthew G. Knepley 
4446858538eSMatthew G. Knepley /*@
445dce8aebaSBarry Smith   DMSetCoordinates - Sets into the `DM` a global vector that holds the coordinates
4466858538eSMatthew G. Knepley 
447356ea6bcSBarry Smith   Logically Collective
4486858538eSMatthew G. Knepley 
4496858538eSMatthew G. Knepley   Input Parameters:
450dce8aebaSBarry Smith + dm - the `DM`
4516858538eSMatthew G. Knepley - c  - coordinate vector
4526858538eSMatthew G. Knepley 
4536858538eSMatthew G. Knepley   Level: intermediate
4546858538eSMatthew G. Knepley 
455dce8aebaSBarry Smith   Notes:
456dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
457dce8aebaSBarry Smith 
45820f4b53cSBarry Smith   The vector `c` can be destroyed after the call
459dce8aebaSBarry Smith 
460dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetCoordinateDM()`, `DMDASetUniformCoordinates()`
4616858538eSMatthew G. Knepley @*/
DMSetCoordinates(DM dm,Vec c)462d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinates(DM dm, Vec c)
463d71ae5a4SJacob Faibussowitsch {
4646858538eSMatthew G. Knepley   PetscFunctionBegin;
4656858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4666858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
4676858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
4686858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
4696858538eSMatthew G. Knepley   dm->coordinates[0].x = c;
4706858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
4716858538eSMatthew G. Knepley   PetscCall(DMCoarsenHookAdd(dm, DMRestrictHook_Coordinates, NULL, NULL));
4726858538eSMatthew G. Knepley   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_Coordinates, NULL, NULL));
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4746858538eSMatthew G. Knepley }
4756858538eSMatthew G. Knepley 
4766858538eSMatthew G. Knepley /*@
477dce8aebaSBarry Smith   DMGetCellCoordinates - Gets a global vector with the cellwise coordinates associated with the `DM`.
4786858538eSMatthew G. Knepley 
47920f4b53cSBarry Smith   Collective
4806858538eSMatthew G. Knepley 
4816858538eSMatthew G. Knepley   Input Parameter:
482dce8aebaSBarry Smith . dm - the `DM`
4836858538eSMatthew G. Knepley 
4846858538eSMatthew G. Knepley   Output Parameter:
4856858538eSMatthew G. Knepley . c - global coordinate vector
4866858538eSMatthew G. Knepley 
487dce8aebaSBarry Smith   Level: intermediate
488dce8aebaSBarry Smith 
489dce8aebaSBarry Smith   Notes:
490dce8aebaSBarry Smith   This is a borrowed reference, so the user should NOT destroy this vector. When the `DM` is
4910b3275a6SBarry Smith   destroyed `c` will no longer be valid.
4926858538eSMatthew G. Knepley 
4936858538eSMatthew G. Knepley   Each process has only the locally-owned portion of the global coordinates (does NOT have the ghost coordinates).
4946858538eSMatthew G. Knepley 
495dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
4966858538eSMatthew G. Knepley @*/
DMGetCellCoordinates(DM dm,Vec * c)497d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinates(DM dm, Vec *c)
498d71ae5a4SJacob Faibussowitsch {
4996858538eSMatthew G. Knepley   PetscFunctionBegin;
5006858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5014f572ea9SToby Isaac   PetscAssertPointer(c, 2);
5026858538eSMatthew G. Knepley   if (!dm->coordinates[1].x && dm->coordinates[1].xl) {
5036858538eSMatthew G. Knepley     DM cdm = NULL;
5046858538eSMatthew G. Knepley 
50511ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
5066858538eSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(cdm, &dm->coordinates[1].x));
5076858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].x, "DG coordinates"));
5086858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalBegin(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
5096858538eSMatthew G. Knepley     PetscCall(DMLocalToGlobalEnd(cdm, dm->coordinates[1].xl, INSERT_VALUES, dm->coordinates[1].x));
5106858538eSMatthew G. Knepley   }
5116858538eSMatthew G. Knepley   *c = dm->coordinates[1].x;
5123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5136858538eSMatthew G. Knepley }
5146858538eSMatthew G. Knepley 
5156858538eSMatthew G. Knepley /*@
516dce8aebaSBarry Smith   DMSetCellCoordinates - Sets into the `DM` a global vector that holds the cellwise coordinates
5176858538eSMatthew G. Knepley 
51820f4b53cSBarry Smith   Collective
5196858538eSMatthew G. Knepley 
5206858538eSMatthew G. Knepley   Input Parameters:
521dce8aebaSBarry Smith + dm - the `DM`
5226858538eSMatthew G. Knepley - c  - cellwise coordinate vector
5236858538eSMatthew G. Knepley 
5246858538eSMatthew G. Knepley   Level: intermediate
5256858538eSMatthew G. Knepley 
526dce8aebaSBarry Smith   Notes:
527dce8aebaSBarry Smith   The coordinates do not include those for ghost points, which are in the local vector.
528dce8aebaSBarry Smith 
52920f4b53cSBarry Smith   The vector `c` should be destroyed by the caller.
530dce8aebaSBarry Smith 
531dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinates()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMGetCellCoordinatesLocal()`, `DMGetCellCoordinateDM()`
5326858538eSMatthew G. Knepley @*/
DMSetCellCoordinates(DM dm,Vec c)533d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinates(DM dm, Vec c)
534d71ae5a4SJacob Faibussowitsch {
5356858538eSMatthew G. Knepley   PetscFunctionBegin;
5366858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5376858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
5386858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
5396858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
5406858538eSMatthew G. Knepley   dm->coordinates[1].x = c;
5416858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5436858538eSMatthew G. Knepley }
5446858538eSMatthew G. Knepley 
5456858538eSMatthew G. Knepley /*@
546dce8aebaSBarry Smith   DMGetCoordinatesLocalSetUp - Prepares a local vector of coordinates, so that `DMGetCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
5476858538eSMatthew G. Knepley 
54820f4b53cSBarry Smith   Collective
5496858538eSMatthew G. Knepley 
5506858538eSMatthew G. Knepley   Input Parameter:
551dce8aebaSBarry Smith . dm - the `DM`
5526858538eSMatthew G. Knepley 
5536858538eSMatthew G. Knepley   Level: advanced
5546858538eSMatthew G. Knepley 
555dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMGetCoordinatesLocalNoncollective()`
5566858538eSMatthew G. Knepley @*/
DMGetCoordinatesLocalSetUp(DM dm)557d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalSetUp(DM dm)
558d71ae5a4SJacob Faibussowitsch {
5596858538eSMatthew G. Knepley   PetscFunctionBegin;
5606858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5616858538eSMatthew G. Knepley   if (!dm->coordinates[0].xl && dm->coordinates[0].x) {
5626858538eSMatthew G. Knepley     DM       cdm = NULL;
5639d92f5ffSMatthew G. Knepley     PetscInt bs;
5646858538eSMatthew G. Knepley 
5656858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
5666858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[0].xl));
56721c42226SMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "Local Coordinates"));
5689d92f5ffSMatthew G. Knepley     // If the size of the vector is 0, it will not get the right block size
5699d92f5ffSMatthew G. Knepley     PetscCall(VecGetBlockSize(dm->coordinates[0].x, &bs));
5709d92f5ffSMatthew G. Knepley     PetscCall(VecSetBlockSize(dm->coordinates[0].xl, bs));
5716858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[0].xl, "coordinates"));
5726858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5736858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[0].x, INSERT_VALUES, dm->coordinates[0].xl));
5746858538eSMatthew G. Knepley   }
5753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5766858538eSMatthew G. Knepley }
5776858538eSMatthew G. Knepley 
5786858538eSMatthew G. Knepley /*@
579dce8aebaSBarry Smith   DMGetCoordinatesLocal - Gets a local vector with the coordinates associated with the `DM`.
5806858538eSMatthew G. Knepley 
58120f4b53cSBarry Smith   Collective the first time it is called
5826858538eSMatthew G. Knepley 
5836858538eSMatthew G. Knepley   Input Parameter:
584dce8aebaSBarry Smith . dm - the `DM`
5856858538eSMatthew G. Knepley 
5866858538eSMatthew G. Knepley   Output Parameter:
5876858538eSMatthew G. Knepley . c - coordinate vector
5886858538eSMatthew G. Knepley 
589dce8aebaSBarry Smith   Level: intermediate
590dce8aebaSBarry Smith 
591dce8aebaSBarry Smith   Notes:
5920b3275a6SBarry Smith   This is a borrowed reference, so the user should NOT destroy `c`
5936858538eSMatthew G. Knepley 
5946858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
5956858538eSMatthew G. Knepley 
596dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
5976858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
5986858538eSMatthew G. Knepley 
599dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`, `DMGetCoordinatesLocalNoncollective()`
6006858538eSMatthew G. Knepley @*/
DMGetCoordinatesLocal(DM dm,Vec * c)601d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocal(DM dm, Vec *c)
602d71ae5a4SJacob Faibussowitsch {
6036858538eSMatthew G. Knepley   PetscFunctionBegin;
6046858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6054f572ea9SToby Isaac   PetscAssertPointer(c, 2);
6066858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
6076858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6096858538eSMatthew G. Knepley }
6106858538eSMatthew G. Knepley 
6116858538eSMatthew G. Knepley /*@
612dce8aebaSBarry Smith   DMGetCoordinatesLocalNoncollective - Non-collective version of `DMGetCoordinatesLocal()`. Fails if global coordinates have been set and `DMGetCoordinatesLocalSetUp()` not called.
6136858538eSMatthew G. Knepley 
61420f4b53cSBarry Smith   Not Collective
6156858538eSMatthew G. Knepley 
6166858538eSMatthew G. Knepley   Input Parameter:
617dce8aebaSBarry Smith . dm - the `DM`
6186858538eSMatthew G. Knepley 
6196858538eSMatthew G. Knepley   Output Parameter:
6206858538eSMatthew G. Knepley . c - coordinate vector
6216858538eSMatthew G. Knepley 
6226858538eSMatthew G. Knepley   Level: advanced
6236858538eSMatthew G. Knepley 
624dce8aebaSBarry Smith   Note:
625dce8aebaSBarry Smith   A previous call to  `DMGetCoordinatesLocal()` or `DMGetCoordinatesLocalSetUp()` ensures that a call to this function will not error.
626dce8aebaSBarry Smith 
627dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinatesLocal()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6286858538eSMatthew G. Knepley @*/
DMGetCoordinatesLocalNoncollective(DM dm,Vec * c)629d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalNoncollective(DM dm, Vec *c)
630d71ae5a4SJacob Faibussowitsch {
6316858538eSMatthew G. Knepley   PetscFunctionBegin;
6326858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6334f572ea9SToby Isaac   PetscAssertPointer(c, 2);
6346858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[0].xl || !dm->coordinates[0].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called");
6356858538eSMatthew G. Knepley   *c = dm->coordinates[0].xl;
6363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6376858538eSMatthew G. Knepley }
6386858538eSMatthew G. Knepley 
6396858538eSMatthew G. Knepley /*@
640dce8aebaSBarry Smith   DMGetCoordinatesLocalTuple - Gets a local vector with the coordinates of specified points and the section describing its layout.
6416858538eSMatthew G. Knepley 
64220f4b53cSBarry Smith   Not Collective
6436858538eSMatthew G. Knepley 
6446858538eSMatthew G. Knepley   Input Parameters:
645dce8aebaSBarry Smith + dm - the `DM`
646dce8aebaSBarry Smith - p  - the `IS` of points whose coordinates will be returned
6476858538eSMatthew G. Knepley 
6486858538eSMatthew G. Knepley   Output Parameters:
6490b3275a6SBarry Smith + pCoordSection - the `PetscSection` describing the layout of pCoord, i.e. each point corresponds to one point in `p`, and DOFs correspond to coordinates
6500b3275a6SBarry Smith - pCoord        - the `Vec` with coordinates of points in `p`
6516858538eSMatthew G. Knepley 
652dce8aebaSBarry Smith   Level: advanced
653dce8aebaSBarry Smith 
654dce8aebaSBarry Smith   Notes:
655dce8aebaSBarry Smith   `DMGetCoordinatesLocalSetUp()` must be called first. This function employs `DMGetCoordinatesLocalNoncollective()` so it is not collective.
6566858538eSMatthew G. Knepley 
6576858538eSMatthew G. Knepley   This creates a new vector, so the user SHOULD destroy this vector
6586858538eSMatthew G. Knepley 
6596858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
6606858538eSMatthew G. Knepley 
661dce8aebaSBarry Smith   For `DMDA`, in two and three dimensions coordinates are interlaced (x_0,y_0,x_1,y_1,...)
6626858538eSMatthew G. Knepley   and (x_0,y_0,z_0,x_1,y_1,z_1...)
6636858538eSMatthew G. Knepley 
664dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMSetCoordinatesLocal()`, `DMGetCoordinatesLocal()`, `DMGetCoordinatesLocalNoncollective()`, `DMGetCoordinatesLocalSetUp()`, `DMGetCoordinates()`, `DMSetCoordinates()`, `DMGetCoordinateDM()`
6656858538eSMatthew G. Knepley @*/
DMGetCoordinatesLocalTuple(DM dm,IS p,PetscSection * pCoordSection,Vec * pCoord)666d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinatesLocalTuple(DM dm, IS p, PetscSection *pCoordSection, Vec *pCoord)
667d71ae5a4SJacob Faibussowitsch {
6686858538eSMatthew G. Knepley   DM                 cdm;
6696858538eSMatthew G. Knepley   PetscSection       cs, newcs;
6706858538eSMatthew G. Knepley   Vec                coords;
6716858538eSMatthew G. Knepley   const PetscScalar *arr;
6726858538eSMatthew G. Knepley   PetscScalar       *newarr = NULL;
6736858538eSMatthew G. Knepley   PetscInt           n;
6746858538eSMatthew G. Knepley 
6756858538eSMatthew G. Knepley   PetscFunctionBegin;
6766858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6776858538eSMatthew G. Knepley   PetscValidHeaderSpecific(p, IS_CLASSID, 2);
6784f572ea9SToby Isaac   if (pCoordSection) PetscAssertPointer(pCoordSection, 3);
6794f572ea9SToby Isaac   if (pCoord) PetscAssertPointer(pCoord, 4);
6806858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
6816858538eSMatthew G. Knepley   PetscCall(DMGetLocalSection(cdm, &cs));
6826858538eSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
6836858538eSMatthew G. Knepley   PetscCheck(coords, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCoordinatesLocalSetUp() has not been called or coordinates not set");
6846858538eSMatthew G. Knepley   PetscCheck(cdm && cs, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM not supported");
6856858538eSMatthew G. Knepley   PetscCall(VecGetArrayRead(coords, &arr));
6866858538eSMatthew G. Knepley   PetscCall(PetscSectionExtractDofsFromArray(cs, MPIU_SCALAR, arr, p, &newcs, pCoord ? ((void **)&newarr) : NULL));
6876858538eSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(coords, &arr));
6886858538eSMatthew G. Knepley   if (pCoord) {
6896858538eSMatthew G. Knepley     PetscCall(PetscSectionGetStorageSize(newcs, &n));
6906858538eSMatthew G. Knepley     /* set array in two steps to mimic PETSC_OWN_POINTER */
6916858538eSMatthew G. Knepley     PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)p), 1, n, NULL, pCoord));
6926858538eSMatthew G. Knepley     PetscCall(VecReplaceArray(*pCoord, newarr));
6936858538eSMatthew G. Knepley   } else {
6946858538eSMatthew G. Knepley     PetscCall(PetscFree(newarr));
6956858538eSMatthew G. Knepley   }
6969371c9d4SSatish Balay   if (pCoordSection) {
6979371c9d4SSatish Balay     *pCoordSection = newcs;
6989371c9d4SSatish Balay   } else PetscCall(PetscSectionDestroy(&newcs));
6993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7006858538eSMatthew G. Knepley }
7016858538eSMatthew G. Knepley 
7026858538eSMatthew G. Knepley /*@
703dce8aebaSBarry Smith   DMSetCoordinatesLocal - Sets into the `DM` a local vector, including ghost points, that holds the coordinates
7046858538eSMatthew G. Knepley 
70520f4b53cSBarry Smith   Not Collective
7066858538eSMatthew G. Knepley 
7076858538eSMatthew G. Knepley   Input Parameters:
708dce8aebaSBarry Smith + dm - the `DM`
7096858538eSMatthew G. Knepley - c  - coordinate vector
7106858538eSMatthew G. Knepley 
711dce8aebaSBarry Smith   Level: intermediate
712dce8aebaSBarry Smith 
7136858538eSMatthew G. Knepley   Notes:
714dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
715dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
7166858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
7176858538eSMatthew G. Knepley 
7180b3275a6SBarry Smith   The vector `c` should be destroyed by the caller.
7196858538eSMatthew G. Knepley 
720dce8aebaSBarry Smith .seealso: `DM`, `DMGetCoordinatesLocal()`, `DMSetCoordinates()`, `DMGetCoordinates()`, `DMGetCoordinateDM()`
7216858538eSMatthew G. Knepley @*/
DMSetCoordinatesLocal(DM dm,Vec c)722d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinatesLocal(DM dm, Vec c)
723d71ae5a4SJacob Faibussowitsch {
7246858538eSMatthew G. Knepley   PetscFunctionBegin;
7256858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7266858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
7276858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
7286858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].xl));
7296858538eSMatthew G. Knepley   dm->coordinates[0].xl = c;
7306858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
7313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7326858538eSMatthew G. Knepley }
7336858538eSMatthew G. Knepley 
7346858538eSMatthew G. Knepley /*@
735dce8aebaSBarry Smith   DMGetCellCoordinatesLocalSetUp - Prepares a local vector of cellwise coordinates, so that `DMGetCellCoordinatesLocalNoncollective()` can be used as non-collective afterwards.
7366858538eSMatthew G. Knepley 
73720f4b53cSBarry Smith   Collective
7386858538eSMatthew G. Knepley 
7396858538eSMatthew G. Knepley   Input Parameter:
740dce8aebaSBarry Smith . dm - the `DM`
7416858538eSMatthew G. Knepley 
7426858538eSMatthew G. Knepley   Level: advanced
7436858538eSMatthew G. Knepley 
744dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalNoncollective()`
7456858538eSMatthew G. Knepley @*/
DMGetCellCoordinatesLocalSetUp(DM dm)746d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalSetUp(DM dm)
747d71ae5a4SJacob Faibussowitsch {
7486858538eSMatthew G. Knepley   PetscFunctionBegin;
7496858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7506858538eSMatthew G. Knepley   if (!dm->coordinates[1].xl && dm->coordinates[1].x) {
7516858538eSMatthew G. Knepley     DM cdm = NULL;
7526858538eSMatthew G. Knepley 
75311ea91a7SMatthew G. Knepley     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
7546858538eSMatthew G. Knepley     PetscCall(DMCreateLocalVector(cdm, &dm->coordinates[1].xl));
7556858538eSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)dm->coordinates[1].xl, "DG coordinates"));
7566858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalBegin(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7576858538eSMatthew G. Knepley     PetscCall(DMGlobalToLocalEnd(cdm, dm->coordinates[1].x, INSERT_VALUES, dm->coordinates[1].xl));
7586858538eSMatthew G. Knepley   }
7593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7606858538eSMatthew G. Knepley }
7616858538eSMatthew G. Knepley 
7626858538eSMatthew G. Knepley /*@
763dce8aebaSBarry Smith   DMGetCellCoordinatesLocal - Gets a local vector with the cellwise coordinates associated with the `DM`.
7646858538eSMatthew G. Knepley 
76520f4b53cSBarry Smith   Collective
7666858538eSMatthew G. Knepley 
7676858538eSMatthew G. Knepley   Input Parameter:
768dce8aebaSBarry Smith . dm - the `DM`
7696858538eSMatthew G. Knepley 
7706858538eSMatthew G. Knepley   Output Parameter:
7716858538eSMatthew G. Knepley . c - coordinate vector
7726858538eSMatthew G. Knepley 
773dce8aebaSBarry Smith   Level: intermediate
774dce8aebaSBarry Smith 
775dce8aebaSBarry Smith   Notes:
7766858538eSMatthew G. Knepley   This is a borrowed reference, so the user should NOT destroy this vector
7776858538eSMatthew G. Knepley 
7786858538eSMatthew G. Knepley   Each process has the local and ghost coordinates
7796858538eSMatthew G. Knepley 
780dce8aebaSBarry Smith .seealso: `DM`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`, `DMGetCellCoordinatesLocalNoncollective()`
7816858538eSMatthew G. Knepley @*/
DMGetCellCoordinatesLocal(DM dm,Vec * c)782d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocal(DM dm, Vec *c)
783d71ae5a4SJacob Faibussowitsch {
7846858538eSMatthew G. Knepley   PetscFunctionBegin;
7856858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7864f572ea9SToby Isaac   PetscAssertPointer(c, 2);
7876858538eSMatthew G. Knepley   PetscCall(DMGetCellCoordinatesLocalSetUp(dm));
7886858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
7893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7906858538eSMatthew G. Knepley }
7916858538eSMatthew G. Knepley 
7926858538eSMatthew G. Knepley /*@
793dce8aebaSBarry Smith   DMGetCellCoordinatesLocalNoncollective - Non-collective version of `DMGetCellCoordinatesLocal()`. Fails if global cellwise coordinates have been set and `DMGetCellCoordinatesLocalSetUp()` not called.
7946858538eSMatthew G. Knepley 
79520f4b53cSBarry Smith   Not Collective
7966858538eSMatthew G. Knepley 
7976858538eSMatthew G. Knepley   Input Parameter:
798dce8aebaSBarry Smith . dm - the `DM`
7996858538eSMatthew G. Knepley 
8006858538eSMatthew G. Knepley   Output Parameter:
8016858538eSMatthew G. Knepley . c - cellwise coordinate vector
8026858538eSMatthew G. Knepley 
8036858538eSMatthew G. Knepley   Level: advanced
8046858538eSMatthew G. Knepley 
805dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocalSetUp()`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinatesLocal()`, `DMGetCellCoordinates()`, `DMSetCellCoordinates()`, `DMGetCellCoordinateDM()`
8066858538eSMatthew G. Knepley @*/
DMGetCellCoordinatesLocalNoncollective(DM dm,Vec * c)807d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCellCoordinatesLocalNoncollective(DM dm, Vec *c)
808d71ae5a4SJacob Faibussowitsch {
8096858538eSMatthew G. Knepley   PetscFunctionBegin;
8106858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8114f572ea9SToby Isaac   PetscAssertPointer(c, 2);
8126858538eSMatthew G. Knepley   PetscCheck(dm->coordinates[1].xl || !dm->coordinates[1].x, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMGetCellCoordinatesLocalSetUp() has not been called");
8136858538eSMatthew G. Knepley   *c = dm->coordinates[1].xl;
8143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8156858538eSMatthew G. Knepley }
8166858538eSMatthew G. Knepley 
8176858538eSMatthew G. Knepley /*@
818dce8aebaSBarry Smith   DMSetCellCoordinatesLocal - Sets into the `DM` a local vector including ghost points that holds the cellwise coordinates
8196858538eSMatthew G. Knepley 
82020f4b53cSBarry Smith   Not Collective
8216858538eSMatthew G. Knepley 
8226858538eSMatthew G. Knepley   Input Parameters:
823dce8aebaSBarry Smith + dm - the `DM`
8246858538eSMatthew G. Knepley - c  - cellwise coordinate vector
8256858538eSMatthew G. Knepley 
826dce8aebaSBarry Smith   Level: intermediate
827dce8aebaSBarry Smith 
8286858538eSMatthew G. Knepley   Notes:
829dce8aebaSBarry Smith   The coordinates of ghost points can be set using `DMSetCoordinates()`
830dce8aebaSBarry Smith   followed by `DMGetCoordinatesLocal()`. This is intended to enable the
8316858538eSMatthew G. Knepley   setting of ghost coordinates outside of the domain.
8326858538eSMatthew G. Knepley 
8330b3275a6SBarry Smith   The vector `c` should be destroyed by the caller.
8346858538eSMatthew G. Knepley 
835dce8aebaSBarry Smith .seealso: `DM`, `DMGetCellCoordinatesLocal()`, `DMSetCellCoordinates()`, `DMGetCellCoordinates()`, `DMGetCellCoordinateDM()`
8366858538eSMatthew G. Knepley @*/
DMSetCellCoordinatesLocal(DM dm,Vec c)837d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCellCoordinatesLocal(DM dm, Vec c)
838d71ae5a4SJacob Faibussowitsch {
8396858538eSMatthew G. Knepley   PetscFunctionBegin;
8406858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8416858538eSMatthew G. Knepley   if (c) PetscValidHeaderSpecific(c, VEC_CLASSID, 2);
8426858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)c));
8436858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].xl));
8446858538eSMatthew G. Knepley   dm->coordinates[1].xl = c;
8456858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[1].x));
8463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8476858538eSMatthew G. Knepley }
8486858538eSMatthew G. Knepley 
DMGetCoordinateField(DM dm,DMField * field)849d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetCoordinateField(DM dm, DMField *field)
850d71ae5a4SJacob Faibussowitsch {
8516858538eSMatthew G. Knepley   PetscFunctionBegin;
8526858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8534f572ea9SToby Isaac   PetscAssertPointer(field, 2);
854213acdd3SPierre Jolivet   if (!dm->coordinates[0].field) PetscTryTypeMethod(dm, createcoordinatefield, &dm->coordinates[0].field);
8556858538eSMatthew G. Knepley   *field = dm->coordinates[0].field;
8563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8576858538eSMatthew G. Knepley }
8586858538eSMatthew G. Knepley 
DMSetCoordinateField(DM dm,DMField field)859d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetCoordinateField(DM dm, DMField field)
860d71ae5a4SJacob Faibussowitsch {
8616858538eSMatthew G. Knepley   PetscFunctionBegin;
8626858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8636858538eSMatthew G. Knepley   if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2);
8646858538eSMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)field));
8656858538eSMatthew G. Knepley   PetscCall(DMFieldDestroy(&dm->coordinates[0].field));
8666858538eSMatthew G. Knepley   dm->coordinates[0].field = field;
8673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8686858538eSMatthew G. Knepley }
8696858538eSMatthew G. Knepley 
DMSetCellCoordinateField(DM dm,DMField field)8704c712d99Sksagiyam PetscErrorCode DMSetCellCoordinateField(DM dm, DMField field)
8714c712d99Sksagiyam {
8724c712d99Sksagiyam   PetscFunctionBegin;
8734c712d99Sksagiyam   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8744c712d99Sksagiyam   if (field) PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 2);
8754c712d99Sksagiyam   PetscCall(PetscObjectReference((PetscObject)field));
8764c712d99Sksagiyam   PetscCall(DMFieldDestroy(&dm->coordinates[1].field));
8774c712d99Sksagiyam   dm->coordinates[1].field = field;
8784c712d99Sksagiyam   PetscFunctionReturn(PETSC_SUCCESS);
8794c712d99Sksagiyam }
8804c712d99Sksagiyam 
DMGetLocalBoundingBox_Coordinates(DM dm,PetscReal lmin[],PetscReal lmax[],PetscInt cs[],PetscInt ce[])8816c6a6b79SMatthew G. Knepley PetscErrorCode DMGetLocalBoundingBox_Coordinates(DM dm, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[])
882d71ae5a4SJacob Faibussowitsch {
8836858538eSMatthew G. Knepley   Vec         coords = NULL;
8846858538eSMatthew G. Knepley   PetscReal   min[3] = {PETSC_MAX_REAL, PETSC_MAX_REAL, PETSC_MAX_REAL};
8856858538eSMatthew G. Knepley   PetscReal   max[3] = {PETSC_MIN_REAL, PETSC_MIN_REAL, PETSC_MIN_REAL};
8866858538eSMatthew G. Knepley   PetscInt    cdim, i, j;
8876c6a6b79SMatthew G. Knepley   PetscMPIInt size;
8886858538eSMatthew G. Knepley 
8896858538eSMatthew G. Knepley   PetscFunctionBegin;
8906c6a6b79SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
8916858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
8926c6a6b79SMatthew G. Knepley   if (size == 1) {
8936c6a6b79SMatthew G. Knepley     const PetscReal *L, *Lstart;
8946c6a6b79SMatthew G. Knepley 
8956c6a6b79SMatthew G. Knepley     PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
8966c6a6b79SMatthew G. Knepley     if (L) {
8976c6a6b79SMatthew G. Knepley       for (PetscInt d = 0; d < cdim; ++d)
8986c6a6b79SMatthew G. Knepley         if (L[d] > 0.0) {
8996c6a6b79SMatthew G. Knepley           min[d] = Lstart[d];
9006c6a6b79SMatthew G. Knepley           max[d] = Lstart[d] + L[d];
9016c6a6b79SMatthew G. Knepley         }
9026c6a6b79SMatthew G. Knepley     }
9036c6a6b79SMatthew G. Knepley   }
904bdf15c88SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(dm, &coords));
905bdf15c88SMatthew G. Knepley   if (coords) {
906bdf15c88SMatthew G. Knepley     const PetscScalar *local_coords;
907bdf15c88SMatthew G. Knepley     PetscInt           N, Ni;
908bdf15c88SMatthew G. Knepley 
909bdf15c88SMatthew G. Knepley     for (j = cdim; j < 3; ++j) {
910bdf15c88SMatthew G. Knepley       min[j] = 0;
911bdf15c88SMatthew G. Knepley       max[j] = 0;
912bdf15c88SMatthew G. Knepley     }
913bdf15c88SMatthew G. Knepley     PetscCall(VecGetArrayRead(coords, &local_coords));
914bdf15c88SMatthew G. Knepley     PetscCall(VecGetLocalSize(coords, &N));
915bdf15c88SMatthew G. Knepley     Ni = N / cdim;
916bdf15c88SMatthew 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]));
920bdf15c88SMatthew G. Knepley       }
921bdf15c88SMatthew G. Knepley     }
922bdf15c88SMatthew G. Knepley     PetscCall(VecRestoreArrayRead(coords, &local_coords));
923bdf15c88SMatthew G. Knepley     PetscCall(DMGetCellCoordinatesLocal(dm, &coords));
9246858538eSMatthew G. Knepley     if (coords) {
9256858538eSMatthew G. Knepley       PetscCall(VecGetArrayRead(coords, &local_coords));
9266858538eSMatthew G. Knepley       PetscCall(VecGetLocalSize(coords, &N));
9276858538eSMatthew G. Knepley       Ni = N / cdim;
9286858538eSMatthew G. Knepley       for (i = 0; i < Ni; ++i) {
929bdf15c88SMatthew G. Knepley         for (j = 0; j < cdim; ++j) {
930bdf15c88SMatthew G. Knepley           min[j] = PetscMin(min[j], PetscRealPart(local_coords[i * cdim + j]));
931bdf15c88SMatthew G. Knepley           max[j] = PetscMax(max[j], PetscRealPart(local_coords[i * cdim + j]));
9326858538eSMatthew G. Knepley         }
9336858538eSMatthew G. Knepley       }
9346858538eSMatthew G. Knepley       PetscCall(VecRestoreArrayRead(coords, &local_coords));
935bdf15c88SMatthew G. Knepley     }
9366858538eSMatthew G. Knepley     if (lmin) PetscCall(PetscArraycpy(lmin, min, cdim));
9376858538eSMatthew G. Knepley     if (lmax) PetscCall(PetscArraycpy(lmax, max, cdim));
9386c6a6b79SMatthew G. Knepley   }
9396c6a6b79SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
9406c6a6b79SMatthew G. Knepley }
9416c6a6b79SMatthew G. Knepley 
9426c6a6b79SMatthew G. Knepley /*@
9436c6a6b79SMatthew G. Knepley   DMGetLocalBoundingBox - Returns the bounding box for the piece of the `DM` on this process.
9446c6a6b79SMatthew G. Knepley 
9456c6a6b79SMatthew G. Knepley   Not Collective
9466c6a6b79SMatthew G. Knepley 
9476c6a6b79SMatthew G. Knepley   Input Parameter:
9486c6a6b79SMatthew G. Knepley . dm - the `DM`
9496c6a6b79SMatthew G. Knepley 
9506c6a6b79SMatthew G. Knepley   Output Parameters:
9516c6a6b79SMatthew G. Knepley + lmin - local minimum coordinates (length coord dim, optional)
9526c6a6b79SMatthew G. Knepley - lmax - local maximum coordinates (length coord dim, optional)
9536c6a6b79SMatthew G. Knepley 
9546c6a6b79SMatthew G. Knepley   Level: beginner
9556c6a6b79SMatthew G. Knepley 
9566c6a6b79SMatthew G. Knepley   Note:
9576c6a6b79SMatthew G. Knepley   If the `DM` is a `DMDA` and has no coordinates, the index bounds are returned instead.
9586c6a6b79SMatthew G. Knepley 
9596c6a6b79SMatthew G. Knepley .seealso: `DM`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMGetBoundingBox()`
9606c6a6b79SMatthew G. Knepley @*/
DMGetLocalBoundingBox(DM dm,PetscReal lmin[],PetscReal lmax[])9616c6a6b79SMatthew G. Knepley PetscErrorCode DMGetLocalBoundingBox(DM dm, PetscReal lmin[], PetscReal lmax[])
9626c6a6b79SMatthew G. Knepley {
9636c6a6b79SMatthew G. Knepley   PetscFunctionBegin;
9646c6a6b79SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9656c6a6b79SMatthew G. Knepley   PetscUseTypeMethod(dm, getlocalboundingbox, lmin, lmax, NULL, NULL);
9663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9676858538eSMatthew G. Knepley }
9686858538eSMatthew G. Knepley 
9696858538eSMatthew G. Knepley /*@
970dce8aebaSBarry Smith   DMGetBoundingBox - Returns the global bounding box for the `DM`.
9716858538eSMatthew G. Knepley 
9726858538eSMatthew G. Knepley   Collective
9736858538eSMatthew G. Knepley 
9746858538eSMatthew G. Knepley   Input Parameter:
975dce8aebaSBarry Smith . dm - the `DM`
9766858538eSMatthew G. Knepley 
9776858538eSMatthew G. Knepley   Output Parameters:
9786858538eSMatthew G. Knepley + gmin - global minimum coordinates (length coord dim, optional)
97935cb6cd3SPierre Jolivet - gmax - global maximum coordinates (length coord dim, optional)
9806858538eSMatthew G. Knepley 
9816858538eSMatthew G. Knepley   Level: beginner
9826858538eSMatthew G. Knepley 
983dce8aebaSBarry Smith .seealso: `DM`, `DMGetLocalBoundingBox()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`
9846858538eSMatthew G. Knepley @*/
DMGetBoundingBox(DM dm,PetscReal gmin[],PetscReal gmax[])985d71ae5a4SJacob Faibussowitsch PetscErrorCode DMGetBoundingBox(DM dm, PetscReal gmin[], PetscReal gmax[])
986d71ae5a4SJacob Faibussowitsch {
9876858538eSMatthew G. Knepley   PetscReal        lmin[3], lmax[3];
9886c6a6b79SMatthew G. Knepley   const PetscReal *L, *Lstart;
9896858538eSMatthew G. Knepley   PetscInt         cdim;
9906858538eSMatthew G. Knepley 
9916858538eSMatthew G. Knepley   PetscFunctionBegin;
9926858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9936858538eSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
9946858538eSMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox(dm, lmin, lmax));
995e91c04dfSPierre Jolivet   if (gmin) PetscCallMPI(MPIU_Allreduce(lmin, gmin, cdim, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)dm)));
996e91c04dfSPierre Jolivet   if (gmax) PetscCallMPI(MPIU_Allreduce(lmax, gmax, cdim, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)dm)));
9976c6a6b79SMatthew G. Knepley   PetscCall(DMGetPeriodicity(dm, NULL, &Lstart, &L));
9986c6a6b79SMatthew G. Knepley   if (L) {
9996c6a6b79SMatthew G. Knepley     for (PetscInt d = 0; d < cdim; ++d)
10006c6a6b79SMatthew G. Knepley       if (L[d] > 0.0) {
10016c6a6b79SMatthew G. Knepley         gmin[d] = Lstart[d];
10026c6a6b79SMatthew G. Knepley         gmax[d] = Lstart[d] + L[d];
10036c6a6b79SMatthew G. Knepley       }
10046c6a6b79SMatthew G. Knepley   }
10053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10066858538eSMatthew G. Knepley }
10076858538eSMatthew G. Knepley 
DMCreateAffineCoordinates_Internal(DM dm,PetscBool localized)10084c712d99Sksagiyam static PetscErrorCode DMCreateAffineCoordinates_Internal(DM dm, PetscBool localized)
100990612307SJed Brown {
1010e44f6aebSMatthew G. Knepley   DM             cdm;
1011e44f6aebSMatthew G. Knepley   PetscFE        feLinear;
1012e44f6aebSMatthew G. Knepley   DMPolytopeType ct;
10131df12153SMatthew G. Knepley   PetscInt       dim, dE, height, cStart, cEnd, gct;
1014e44f6aebSMatthew G. Knepley 
1015e44f6aebSMatthew G. Knepley   PetscFunctionBegin;
10164c712d99Sksagiyam   if (!localized) {
1017e44f6aebSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdm));
10184c712d99Sksagiyam   } else {
10194c712d99Sksagiyam     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
10204c712d99Sksagiyam   }
10214c712d99Sksagiyam   PetscCheck(cdm, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "No coordinateDM defined");
1022e44f6aebSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
1023e44f6aebSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &dE));
10241df12153SMatthew G. Knepley   PetscCall(DMPlexGetVTKCellHeight(dm, &height));
10251df12153SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, height, &cStart, &cEnd));
1026e44f6aebSMatthew G. Knepley   if (cEnd > cStart) PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1027e44f6aebSMatthew G. Knepley   else ct = DM_POLYTOPE_UNKNOWN;
1028e44f6aebSMatthew G. Knepley   gct = (PetscInt)ct;
1029462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &gct, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm)));
1030e44f6aebSMatthew G. Knepley   ct = (DMPolytopeType)gct;
1031e44f6aebSMatthew G. Knepley   // Work around current bug in PetscDualSpaceSetUp_Lagrange()
1032e44f6aebSMatthew G. Knepley   //   Can be seen in plex_tutorials-ex10_1
1033e44f6aebSMatthew G. Knepley   if (ct != DM_POLYTOPE_SEG_PRISM_TENSOR && ct != DM_POLYTOPE_TRI_PRISM_TENSOR && ct != DM_POLYTOPE_QUAD_PRISM_TENSOR) {
1034e44f6aebSMatthew G. Knepley     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, dE, ct, 1, -1, &feLinear));
10354c712d99Sksagiyam     if (localized) {
10364c712d99Sksagiyam       PetscFE dgfe = NULL;
10374c712d99Sksagiyam 
10384c712d99Sksagiyam       PetscCall(PetscFECreateBrokenElement(feLinear, &dgfe));
10394c712d99Sksagiyam       PetscCall(PetscFEDestroy(&feLinear));
10404c712d99Sksagiyam       feLinear = dgfe;
10414c712d99Sksagiyam     }
1042e44f6aebSMatthew G. Knepley     PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)feLinear));
1043e44f6aebSMatthew G. Knepley     PetscCall(PetscFEDestroy(&feLinear));
1044e44f6aebSMatthew G. Knepley     PetscCall(DMCreateDS(cdm));
1045e44f6aebSMatthew G. Knepley   }
1046e44f6aebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1047e44f6aebSMatthew G. Knepley }
1048e44f6aebSMatthew G. Knepley 
DMGetCoordinateDegree_Internal(DM dm,PetscInt * degree)1049e44f6aebSMatthew G. Knepley PetscErrorCode DMGetCoordinateDegree_Internal(DM dm, PetscInt *degree)
1050e44f6aebSMatthew G. Knepley {
1051e44f6aebSMatthew G. Knepley   DM           cdm;
1052e44f6aebSMatthew G. Knepley   PetscFE      fe;
1053e44f6aebSMatthew G. Knepley   PetscSpace   sp;
1054e44f6aebSMatthew G. Knepley   PetscClassId id;
1055e44f6aebSMatthew G. Knepley 
1056e44f6aebSMatthew G. Knepley   PetscFunctionBegin;
1057e44f6aebSMatthew G. Knepley   *degree = 1;
1058e44f6aebSMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &cdm));
1059e44f6aebSMatthew G. Knepley   PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe));
1060e44f6aebSMatthew G. Knepley   PetscCall(PetscObjectGetClassId((PetscObject)fe, &id));
1061e44f6aebSMatthew G. Knepley   if (id != PETSCFE_CLASSID) PetscFunctionReturn(PETSC_SUCCESS);
1062e44f6aebSMatthew G. Knepley   PetscCall(PetscFEGetBasisSpace(fe, &sp));
1063e44f6aebSMatthew G. Knepley   PetscCall(PetscSpaceGetDegree(sp, degree, NULL));
1064e44f6aebSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
106590612307SJed Brown }
106690612307SJed Brown 
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[])106790df3356SJames Wright static void evaluate_coordinates(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar xnew[])
106890df3356SJames Wright {
106990df3356SJames Wright   for (PetscInt i = 0; i < dim; i++) xnew[i] = x[i];
107090df3356SJames Wright }
107190df3356SJames Wright 
10726858538eSMatthew G. Knepley /*@
1073e44f6aebSMatthew G. Knepley   DMSetCoordinateDisc - Set a coordinate space
10746858538eSMatthew G. Knepley 
10756858538eSMatthew G. Knepley   Input Parameters:
1076dce8aebaSBarry Smith + dm        - The `DM` object
10774688473bSSatish Balay . disc      - The new coordinate discretization or `NULL` to ensure a coordinate discretization exists
10784c712d99Sksagiyam . localized - Set a localized (DG) coordinate space
1079e44f6aebSMatthew G. Knepley - project   - Project coordinates to new discretization
10806858538eSMatthew G. Knepley 
10816858538eSMatthew G. Knepley   Level: intermediate
10826858538eSMatthew G. Knepley 
1083dce8aebaSBarry Smith   Notes:
1084e44f6aebSMatthew G. Knepley   A `PetscFE` defines an approximation space using a `PetscSpace`, which represents the basis functions, and a `PetscDualSpace`, which defines the interpolation operation in the space.
1085dce8aebaSBarry Smith 
1086dce8aebaSBarry Smith   This function takes the current mesh coordinates, which are discretized using some `PetscFE` space, and projects this function into a new `PetscFE` space.
1087e44f6aebSMatthew G. Knepley   The coordinate projection is done on the continuous coordinates, but the discontinuous coordinates are not updated.
1088dce8aebaSBarry Smith 
10890b3275a6SBarry Smith   Developer Note:
1090dce8aebaSBarry Smith   With more effort, we could directly project the discontinuous coordinates also.
1091dce8aebaSBarry Smith 
1092dce8aebaSBarry Smith .seealso: `DM`, `PetscFE`, `DMGetCoordinateField()`
10936858538eSMatthew G. Knepley @*/
DMSetCoordinateDisc(DM dm,PetscFE disc,PetscBool localized,PetscBool project)10944c712d99Sksagiyam PetscErrorCode DMSetCoordinateDisc(DM dm, PetscFE disc, PetscBool localized, PetscBool project)
1095d71ae5a4SJacob Faibussowitsch {
1096e44f6aebSMatthew G. Knepley   DM           cdmOld, cdmNew;
10976858538eSMatthew G. Knepley   PetscFE      discOld;
10986858538eSMatthew G. Knepley   PetscClassId classid;
10996858538eSMatthew G. Knepley   PetscBool    same_space = PETSC_TRUE;
1100dd4c3f67SMatthew G. Knepley   const char  *prefix;
11016858538eSMatthew G. Knepley 
11026858538eSMatthew G. Knepley   PetscFunctionBegin;
11036858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11046858538eSMatthew G. Knepley   if (disc) PetscValidHeaderSpecific(disc, PETSCFE_CLASSID, 2);
11056858538eSMatthew G. Knepley 
11064c712d99Sksagiyam   /* Note that plexgmsh.c can pass DG element with localized = PETSC_FALSE. */
11074c712d99Sksagiyam   if (!localized) {
11086858538eSMatthew G. Knepley     PetscCall(DMGetCoordinateDM(dm, &cdmOld));
11094c712d99Sksagiyam   } else {
11104c712d99Sksagiyam     PetscCall(DMGetCellCoordinateDM(dm, &cdmOld));
11114c712d99Sksagiyam     if (!cdmOld) {
11124c712d99Sksagiyam       PetscUseTypeMethod(dm, createcellcoordinatedm, &cdmOld);
11134c712d99Sksagiyam       PetscCall(DMSetCellCoordinateDM(dm, cdmOld));
11144c712d99Sksagiyam       PetscCall(DMDestroy(&cdmOld));
11154c712d99Sksagiyam       PetscCall(DMGetCellCoordinateDM(dm, &cdmOld));
11164c712d99Sksagiyam     }
11174c712d99Sksagiyam   }
11186858538eSMatthew G. Knepley   /* Check current discretization is compatible */
11196858538eSMatthew G. Knepley   PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
11206858538eSMatthew G. Knepley   PetscCall(PetscObjectGetClassId((PetscObject)discOld, &classid));
11216858538eSMatthew G. Knepley   if (classid != PETSCFE_CLASSID) {
11226858538eSMatthew G. Knepley     if (classid == PETSC_CONTAINER_CLASSID) {
11234c712d99Sksagiyam       PetscCall(DMCreateAffineCoordinates_Internal(dm, localized));
11246858538eSMatthew G. Knepley       PetscCall(DMGetField(cdmOld, 0, NULL, (PetscObject *)&discOld));
11256858538eSMatthew G. Knepley     } else {
11266858538eSMatthew G. Knepley       const char *discname;
11276858538eSMatthew G. Knepley 
11286858538eSMatthew G. Knepley       PetscCall(PetscObjectGetType((PetscObject)discOld, &discname));
11296858538eSMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)discOld), PETSC_ERR_SUP, "Discretization type %s not supported", discname);
11306858538eSMatthew G. Knepley     }
11316858538eSMatthew G. Knepley   }
1132e44f6aebSMatthew G. Knepley   // Linear space has been created by now
11333ba16761SJacob Faibussowitsch   if (!disc) PetscFunctionReturn(PETSC_SUCCESS);
1134e44f6aebSMatthew G. Knepley   // Check if the new space is the same as the old modulo quadrature
1135e44f6aebSMatthew G. Knepley   {
11366858538eSMatthew G. Knepley     PetscDualSpace dsOld, ds;
11376858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(discOld, &dsOld));
11386858538eSMatthew G. Knepley     PetscCall(PetscFEGetDualSpace(disc, &ds));
11396858538eSMatthew G. Knepley     PetscCall(PetscDualSpaceEqual(dsOld, ds, &same_space));
11406858538eSMatthew G. Knepley   }
1141e44f6aebSMatthew G. Knepley   // Make a fresh clone of the coordinate DM
11426858538eSMatthew G. Knepley   PetscCall(DMClone(cdmOld, &cdmNew));
1143dd4c3f67SMatthew G. Knepley   cdmNew->cloneOpts = PETSC_TRUE;
1144dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdmOld, &prefix));
1145dd4c3f67SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)cdmNew, prefix));
11466858538eSMatthew G. Knepley   PetscCall(DMSetField(cdmNew, 0, NULL, (PetscObject)disc));
11476858538eSMatthew G. Knepley   PetscCall(DMCreateDS(cdmNew));
1148e44f6aebSMatthew G. Knepley   {
1149e44f6aebSMatthew G. Knepley     PetscDS ds, nds;
1150e44f6aebSMatthew G. Knepley 
1151e44f6aebSMatthew G. Knepley     PetscCall(DMGetDS(cdmOld, &ds));
1152e44f6aebSMatthew G. Knepley     PetscCall(DMGetDS(cdmNew, &nds));
1153e44f6aebSMatthew G. Knepley     PetscCall(PetscDSCopyConstants(ds, nds));
1154e44f6aebSMatthew G. Knepley   }
115581316d9bSJames Wright   if (cdmOld->periodic.setup) {
115681316d9bSJames Wright     PetscSF dummy;
115781316d9bSJames Wright     // Force IsoperiodicPointSF to be built, required for periodic coordinate setup
115881316d9bSJames Wright     PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &dummy));
115981316d9bSJames Wright     cdmNew->periodic.setup = cdmOld->periodic.setup;
116081316d9bSJames Wright     PetscCall(cdmNew->periodic.setup(cdmNew));
116181316d9bSJames Wright   }
1162dd4c3f67SMatthew G. Knepley   if (dm->setfromoptionscalled) PetscCall(DMSetFromOptions(cdmNew));
1163e44f6aebSMatthew G. Knepley   if (project) {
116453b735f8SJames Wright     Vec      coordsOld, coordsNew;
116581316d9bSJames Wright     PetscInt num_face_sfs = 0;
116653b735f8SJames Wright 
116781316d9bSJames Wright     PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, NULL));
116881316d9bSJames Wright     if (num_face_sfs) { // Isoperiodicity requires projecting the local coordinates
116953b735f8SJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &coordsOld));
117053b735f8SJames Wright       PetscCall(DMCreateLocalVector(cdmNew, &coordsNew));
117153b735f8SJames Wright       PetscCall(PetscObjectSetName((PetscObject)coordsNew, "coordinates"));
1172f0ac6e35SMatthew G. Knepley       if (same_space) {
1173f0ac6e35SMatthew G. Knepley         // Need to copy so that the new vector has the right dm
1174f0ac6e35SMatthew G. Knepley         PetscCall(VecCopy(coordsOld, coordsNew));
1175e44f6aebSMatthew G. Knepley       } else {
117690df3356SJames Wright         void (*funcs[])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {evaluate_coordinates};
117790df3356SJames Wright 
117890df3356SJames Wright         // We can't call DMProjectField directly because it depends on KSP for DMGlobalToLocalSolve(), but we can use the core strategy
117990df3356SJames Wright         PetscCall(DMSetCoordinateDM(cdmNew, cdmOld));
118090df3356SJames Wright         // See DMPlexRemapGeometry() for a similar pattern handling the coordinate field
118190df3356SJames Wright         DMField cf;
118290df3356SJames Wright         PetscCall(DMGetCoordinateField(dm, &cf));
118390df3356SJames Wright         cdmNew->coordinates[0].field = cf;
118453b735f8SJames Wright         PetscCall(DMProjectFieldLocal(cdmNew, 0.0, NULL, funcs, INSERT_VALUES, coordsNew));
118590df3356SJames Wright         cdmNew->coordinates[0].field = NULL;
118690df3356SJames Wright         PetscCall(DMSetCoordinateDM(cdmNew, NULL));
118753b735f8SJames Wright       }
118853b735f8SJames Wright       PetscCall(DMSetCoordinatesLocal(dm, coordsNew));
118953b735f8SJames Wright       PetscCall(VecDestroy(&coordsNew));
119053b735f8SJames Wright     } else {
119153b735f8SJames Wright       PetscCall(DMGetCoordinates(dm, &coordsOld));
119253b735f8SJames Wright       PetscCall(DMCreateGlobalVector(cdmNew, &coordsNew));
119353b735f8SJames Wright       if (same_space) {
119453b735f8SJames Wright         // Need to copy so that the new vector has the right dm
119553b735f8SJames Wright         PetscCall(VecCopy(coordsOld, coordsNew));
119690df3356SJames Wright       } else {
1197e44f6aebSMatthew G. Knepley         Mat In;
1198e44f6aebSMatthew G. Knepley 
1199e44f6aebSMatthew G. Knepley         PetscCall(DMCreateInterpolation(cdmOld, cdmNew, &In, NULL));
1200e44f6aebSMatthew G. Knepley         PetscCall(MatMult(In, coordsOld, coordsNew));
1201e44f6aebSMatthew G. Knepley         PetscCall(MatDestroy(&In));
1202e44f6aebSMatthew G. Knepley       }
1203e44f6aebSMatthew G. Knepley       PetscCall(DMSetCoordinates(dm, coordsNew));
1204e44f6aebSMatthew G. Knepley       PetscCall(VecDestroy(&coordsNew));
12056858538eSMatthew G. Knepley     }
120653b735f8SJames Wright   }
12076858538eSMatthew G. Knepley   /* Set new coordinate structures */
12084c712d99Sksagiyam   if (!localized) {
12096858538eSMatthew G. Knepley     PetscCall(DMSetCoordinateField(dm, NULL));
12106858538eSMatthew G. Knepley     PetscCall(DMSetCoordinateDM(dm, cdmNew));
12114c712d99Sksagiyam   } else {
12124c712d99Sksagiyam     PetscCall(DMSetCellCoordinateField(dm, NULL));
12134c712d99Sksagiyam     PetscCall(DMSetCellCoordinateDM(dm, cdmNew));
12144c712d99Sksagiyam   }
12156858538eSMatthew G. Knepley   PetscCall(DMDestroy(&cdmNew));
12163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12176858538eSMatthew G. Knepley }
12186858538eSMatthew G. Knepley 
12196858538eSMatthew G. Knepley /*@
12200b3275a6SBarry Smith   DMLocatePoints - Locate the points in `v` in the mesh and return a `PetscSF` of the containing cells
12216858538eSMatthew G. Knepley 
122220f4b53cSBarry Smith   Collective
12236858538eSMatthew G. Knepley 
12246858538eSMatthew G. Knepley   Input Parameters:
1225dce8aebaSBarry Smith + dm    - The `DM`
1226dce8aebaSBarry Smith - ltype - The type of point location, e.g. `DM_POINTLOCATION_NONE` or `DM_POINTLOCATION_NEAREST`
12276858538eSMatthew G. Knepley 
12286858538eSMatthew G. Knepley   Input/Output Parameters:
1229dce8aebaSBarry Smith + v      - The `Vec` of points, on output contains the nearest mesh points to the given points if `DM_POINTLOCATION_NEAREST` is used
123020f4b53cSBarry Smith - cellSF - Points to either `NULL`, or a `PetscSF` with guesses for which cells contain each point;
12310b3275a6SBarry Smith            on output, the `PetscSF` containing the MPI ranks and local indices of the containing points
12326858538eSMatthew G. Knepley 
12336858538eSMatthew G. Knepley   Level: developer
12346858538eSMatthew G. Knepley 
12356858538eSMatthew G. Knepley   Notes:
12360b3275a6SBarry Smith   To do a search of the local cells of the mesh, `v` should have `PETSC_COMM_SELF` as its communicator.
12370b3275a6SBarry Smith   To do a search of all the cells in the distributed mesh, `v` should have the same MPI communicator as `dm`.
12386858538eSMatthew G. Knepley 
1239d8206211SMatthew G. Knepley   Points will only be located in owned cells, not overlap cells arising from `DMPlexDistribute()` or other overlapping distributions.
1240d8206211SMatthew G. Knepley 
124120f4b53cSBarry Smith   If *cellSF is `NULL` on input, a `PetscSF` will be created.
124220f4b53cSBarry Smith   If *cellSF is not `NULL` on input, it should point to an existing `PetscSF`, whose graph will be used as initial guesses.
12436858538eSMatthew G. Knepley 
12446858538eSMatthew G. Knepley   An array that maps each point to its containing cell can be obtained with
1245dce8aebaSBarry Smith .vb
1246dce8aebaSBarry Smith     const PetscSFNode *cells;
1247dce8aebaSBarry Smith     PetscInt           nFound;
1248dce8aebaSBarry Smith     const PetscInt    *found;
12496858538eSMatthew G. Knepley 
1250dce8aebaSBarry Smith     PetscSFGetGraph(cellSF,NULL,&nFound,&found,&cells);
1251dce8aebaSBarry Smith .ve
12526858538eSMatthew G. Knepley 
12530b3275a6SBarry Smith   Where cells[i].rank is the MPI rank of the process owning the cell containing point found[i] (or i if found == NULL), and cells[i].index is
12540b3275a6SBarry Smith   the index of the cell in its MPI process' local numbering. This rank is in the communicator for `v`, so if `v` is on `PETSC_COMM_SELF` then the rank will always be 0.
12556858538eSMatthew G. Knepley 
1256dce8aebaSBarry Smith .seealso: `DM`, `DMSetCoordinates()`, `DMSetCoordinatesLocal()`, `DMGetCoordinates()`, `DMGetCoordinatesLocal()`, `DMPointLocationType`
12576858538eSMatthew G. Knepley @*/
DMLocatePoints(DM dm,Vec v,DMPointLocationType ltype,PetscSF * cellSF)1258d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints(DM dm, Vec v, DMPointLocationType ltype, PetscSF *cellSF)
1259d71ae5a4SJacob Faibussowitsch {
12606858538eSMatthew G. Knepley   PetscFunctionBegin;
12616858538eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12626858538eSMatthew G. Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 2);
12634f572ea9SToby Isaac   PetscAssertPointer(cellSF, 4);
12646858538eSMatthew G. Knepley   if (*cellSF) {
12656858538eSMatthew G. Knepley     PetscMPIInt result;
12666858538eSMatthew G. Knepley 
12676858538eSMatthew G. Knepley     PetscValidHeaderSpecific(*cellSF, PETSCSF_CLASSID, 4);
12686858538eSMatthew G. Knepley     PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)v), PetscObjectComm((PetscObject)*cellSF), &result));
12696858538eSMatthew 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");
12706858538eSMatthew G. Knepley   } else {
12716858538eSMatthew G. Knepley     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)v), cellSF));
12726858538eSMatthew G. Knepley   }
12736858538eSMatthew G. Knepley   PetscCall(PetscLogEventBegin(DM_LocatePoints, dm, 0, 0, 0));
1274dbbe0bcdSBarry Smith   PetscUseTypeMethod(dm, locatepoints, v, ltype, *cellSF);
12756858538eSMatthew G. Knepley   PetscCall(PetscLogEventEnd(DM_LocatePoints, dm, 0, 0, 0));
12763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12776858538eSMatthew G. Knepley }
1278