xref: /petsc/src/dm/impls/plex/tests/ex17.c (revision 274abdb237baef925d63e54f8ba8d7d03a350e6e)
1c4762a1bSJed Brown static char help[] = "Tests for point location\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsf.h>
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown 
6*274abdb2SMatthew G. Knepley typedef struct {
7*274abdb2SMatthew G. Knepley   PetscBool centroids;
8*274abdb2SMatthew G. Knepley   PetscBool custom;
9*274abdb2SMatthew G. Knepley } AppCtx;
10*274abdb2SMatthew G. Knepley 
11*274abdb2SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12*274abdb2SMatthew G. Knepley {
13*274abdb2SMatthew G. Knepley   PetscFunctionBeginUser;
14*274abdb2SMatthew G. Knepley   options->centroids = PETSC_TRUE;
15*274abdb2SMatthew G. Knepley   options->custom    = PETSC_FALSE;
16*274abdb2SMatthew G. Knepley 
17*274abdb2SMatthew G. Knepley   PetscOptionsBegin(comm, "", "Point Location Options", "DMPLEX");
18*274abdb2SMatthew G. Knepley   PetscCall(PetscOptionsBool("-centroids", "Locate cell centroids", "ex17.c", options->centroids, &options->centroids, NULL));
19*274abdb2SMatthew G. Knepley   PetscCall(PetscOptionsBool("-custom", "Locate user-defined points", "ex17.c", options->custom, &options->custom, NULL));
20*274abdb2SMatthew G. Knepley   PetscOptionsEnd();
21*274abdb2SMatthew G. Knepley   PetscFunctionReturn(0);
22*274abdb2SMatthew G. Knepley }
23*274abdb2SMatthew G. Knepley 
24d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown   PetscFunctionBeginUser;
279566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
289566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
299566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
309566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
31c4762a1bSJed Brown   PetscFunctionReturn(0);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34*274abdb2SMatthew G. Knepley static PetscErrorCode TestCentroidLocation(DM dm, AppCtx *user)
35d71ae5a4SJacob Faibussowitsch {
363285882aSMatthew G. Knepley   Vec                points;
373285882aSMatthew G. Knepley   PetscSF            cellSF = NULL;
383285882aSMatthew G. Knepley   const PetscSFNode *cells;
393285882aSMatthew G. Knepley   PetscScalar       *a;
403285882aSMatthew G. Knepley   PetscInt           cdim, n;
41c4762a1bSJed Brown   PetscInt           cStart, cEnd, c;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   PetscFunctionBeginUser;
44*274abdb2SMatthew G. Knepley   if (!user->centroids) PetscFunctionReturn(0);
459566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
468fb5bd83SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
48c4762a1bSJed Brown   /* Locate all centroids */
499566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart) * cdim, &points));
509566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(points, cdim));
519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(points, &a));
52c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
53c4762a1bSJed Brown     PetscReal centroid[3];
543285882aSMatthew G. Knepley     PetscInt  off = (c - cStart) * cdim, d;
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
573285882aSMatthew G. Knepley     for (d = 0; d < cdim; ++d) a[off + d] = centroid[d];
58c4762a1bSJed Brown   }
599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(points, &a));
609566063dSJacob Faibussowitsch   PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF));
619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&points));
629566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(cellSF, NULL, &n, NULL, &cells));
633285882aSMatthew G. Knepley   if (n != (cEnd - cStart)) {
643285882aSMatthew G. Knepley     for (c = 0; c < n; ++c) {
6563a3b9bcSJacob Faibussowitsch       if (cells[c].index != c + cStart) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Could not locate centroid of cell %" PetscInt_FMT ", error %" PetscInt_FMT "\n", c + cStart, cells[c].index));
663285882aSMatthew G. Knepley     }
6763a3b9bcSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Located %" PetscInt_FMT " points instead of %" PetscInt_FMT, n, cEnd - cStart);
683285882aSMatthew G. Knepley   }
69ad540459SPierre Jolivet   for (c = cStart; c < cEnd; ++c) PetscCheck(cells[c - cStart].index == c, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not locate centroid of cell %" PetscInt_FMT ", instead found %" PetscInt_FMT, c, cells[c - cStart].index);
709566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&cellSF));
71c4762a1bSJed Brown   PetscFunctionReturn(0);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74*274abdb2SMatthew G. Knepley static PetscErrorCode TestCustomLocation(DM dm, AppCtx *user)
75*274abdb2SMatthew G. Knepley {
76*274abdb2SMatthew G. Knepley   PetscSF            cellSF = NULL;
77*274abdb2SMatthew G. Knepley   const PetscSFNode *cells;
78*274abdb2SMatthew G. Knepley   const PetscInt    *found;
79*274abdb2SMatthew G. Knepley   Vec                points;
80*274abdb2SMatthew G. Knepley   PetscScalar        coords[2] = {0.5, 0.5};
81*274abdb2SMatthew G. Knepley   PetscInt           cdim, Np = 1, Nfd;
82*274abdb2SMatthew G. Knepley   PetscMPIInt        rank;
83*274abdb2SMatthew G. Knepley   MPI_Comm           comm;
84*274abdb2SMatthew G. Knepley 
85*274abdb2SMatthew G. Knepley   PetscFunctionBeginUser;
86*274abdb2SMatthew G. Knepley   if (!user->custom) PetscFunctionReturn(0);
87*274abdb2SMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
88*274abdb2SMatthew G. Knepley 
89*274abdb2SMatthew G. Knepley   // Locate serially on each process
90*274abdb2SMatthew G. Knepley   PetscCall(VecCreate(PETSC_COMM_SELF, &points));
91*274abdb2SMatthew G. Knepley   PetscCall(VecSetBlockSize(points, cdim));
92*274abdb2SMatthew G. Knepley   PetscCall(VecSetSizes(points, Np * cdim, PETSC_DETERMINE));
93*274abdb2SMatthew G. Knepley   PetscCall(VecSetFromOptions(points));
94*274abdb2SMatthew G. Knepley   for (PetscInt p = 0; p < Np; ++p) {
95*274abdb2SMatthew G. Knepley     const PetscInt idx[2] = {p * cdim, p * cdim + 1};
96*274abdb2SMatthew G. Knepley     PetscCall(VecSetValues(points, cdim, idx, coords, INSERT_VALUES));
97*274abdb2SMatthew G. Knepley   }
98*274abdb2SMatthew G. Knepley   PetscCall(VecAssemblyBegin(points));
99*274abdb2SMatthew G. Knepley   PetscCall(VecAssemblyEnd(points));
100*274abdb2SMatthew G. Knepley 
101*274abdb2SMatthew G. Knepley   PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF));
102*274abdb2SMatthew G. Knepley 
103*274abdb2SMatthew G. Knepley   PetscCall(PetscSFGetGraph(cellSF, NULL, &Nfd, &found, &cells));
104*274abdb2SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
105*274abdb2SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(comm, &rank));
106*274abdb2SMatthew G. Knepley   PetscCall(PetscSynchronizedPrintf(comm, "[%d] Found %" PetscInt_FMT " particles\n", rank, Nfd));
107*274abdb2SMatthew G. Knepley   for (PetscInt p = 0; p < Nfd; ++p) {
108*274abdb2SMatthew G. Knepley     const PetscInt     point = found ? found[p] : p;
109*274abdb2SMatthew G. Knepley     const PetscScalar *array;
110*274abdb2SMatthew G. Knepley     PetscScalar       *ccoords = NULL;
111*274abdb2SMatthew G. Knepley     PetscInt           numCoords;
112*274abdb2SMatthew G. Knepley     PetscBool          isDG;
113*274abdb2SMatthew G. Knepley 
114*274abdb2SMatthew G. Knepley     // Since the v comm is SELF, rank is always 0
115*274abdb2SMatthew G. Knepley     PetscCall(PetscSynchronizedPrintf(comm, "  point %" PetscInt_FMT " cell %" PetscInt_FMT "\n", point, cells[p].index));
116*274abdb2SMatthew G. Knepley     PetscCall(DMPlexGetCellCoordinates(dm, cells[p].index, &isDG, &numCoords, &array, &ccoords));
117*274abdb2SMatthew G. Knepley     for (PetscInt c = 0; c < numCoords / cdim; ++c) {
118*274abdb2SMatthew G. Knepley       PetscCall(PetscSynchronizedPrintf(comm, "  "));
119*274abdb2SMatthew G. Knepley       for (PetscInt d = 0; d < cdim; ++d) PetscCall(PetscSynchronizedPrintf(comm, " %g", (double)PetscRealPart(ccoords[c * cdim + d])));
120*274abdb2SMatthew G. Knepley       PetscCall(PetscSynchronizedPrintf(comm, "\n"));
121*274abdb2SMatthew G. Knepley     }
122*274abdb2SMatthew G. Knepley     PetscCall(DMPlexRestoreCellCoordinates(dm, cells[p].index, &isDG, &numCoords, &array, &ccoords));
123*274abdb2SMatthew G. Knepley   }
124*274abdb2SMatthew G. Knepley   PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
125*274abdb2SMatthew G. Knepley 
126*274abdb2SMatthew G. Knepley   PetscCall(PetscSFDestroy(&cellSF));
127*274abdb2SMatthew G. Knepley   PetscCall(VecDestroy(&points));
128*274abdb2SMatthew G. Knepley   PetscFunctionReturn(0);
129*274abdb2SMatthew G. Knepley }
130*274abdb2SMatthew G. Knepley 
131d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
132d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown   DM     dm;
134*274abdb2SMatthew G. Knepley   AppCtx user;
135c4762a1bSJed Brown 
136327415f7SBarry Smith   PetscFunctionBeginUser;
1379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
138*274abdb2SMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1399566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
140*274abdb2SMatthew G. Knepley   PetscCall(TestCentroidLocation(dm, &user));
141*274abdb2SMatthew G. Knepley   PetscCall(TestCustomLocation(dm, &user));
1429566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
144b122ec5aSJacob Faibussowitsch   return 0;
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown /*TEST
148c4762a1bSJed Brown 
149b26b5bf9SMatthew G. Knepley   testset:
150b26b5bf9SMatthew G. Knepley     args: -dm_plex_dim 1 -dm_plex_box_faces 10
151b26b5bf9SMatthew G. Knepley 
152c4762a1bSJed Brown     test:
1531af33867SMatthew G. Knepley       suffix: seg
154b26b5bf9SMatthew G. Knepley 
155b26b5bf9SMatthew G. Knepley     test:
156b26b5bf9SMatthew G. Knepley       suffix: seg_hash
157b26b5bf9SMatthew G. Knepley       args: -dm_refine 2 -dm_plex_hash_location
1581af33867SMatthew G. Knepley 
1593285882aSMatthew G. Knepley   testset:
1603285882aSMatthew G. Knepley     args: -dm_plex_box_faces 5,5
1613285882aSMatthew G. Knepley 
1621af33867SMatthew G. Knepley     test:
1631af33867SMatthew G. Knepley       suffix: tri
164c4762a1bSJed Brown       requires: triangle
1653285882aSMatthew G. Knepley 
1663285882aSMatthew G. Knepley     test:
1673285882aSMatthew G. Knepley       suffix: tri_hash
1683285882aSMatthew G. Knepley       requires: triangle
1693285882aSMatthew G. Knepley       args: -dm_refine 2 -dm_plex_hash_location
1701af33867SMatthew G. Knepley 
1711af33867SMatthew G. Knepley     test:
1721af33867SMatthew G. Knepley       suffix: quad
1733285882aSMatthew G. Knepley       args: -dm_plex_simplex 0
1743285882aSMatthew G. Knepley 
1753285882aSMatthew G. Knepley     test:
1763285882aSMatthew G. Knepley       suffix: quad_hash
1773285882aSMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_refine 2 -dm_plex_hash_location
1783285882aSMatthew G. Knepley 
1793285882aSMatthew G. Knepley   testset:
1803285882aSMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3
1811af33867SMatthew G. Knepley 
1821af33867SMatthew G. Knepley     test:
1831af33867SMatthew G. Knepley       suffix: tet
1841af33867SMatthew G. Knepley       requires: ctetgen
1853285882aSMatthew G. Knepley 
1863285882aSMatthew G. Knepley     test:
1873285882aSMatthew G. Knepley       suffix: tet_hash
1883285882aSMatthew G. Knepley       requires: ctetgen
1893285882aSMatthew G. Knepley       args: -dm_refine 1 -dm_plex_hash_location
1901af33867SMatthew G. Knepley 
1911af33867SMatthew G. Knepley     test:
1921af33867SMatthew G. Knepley       suffix: hex
1933285882aSMatthew G. Knepley       args: -dm_plex_simplex 0
1943285882aSMatthew G. Knepley 
1953285882aSMatthew G. Knepley     test:
1963285882aSMatthew G. Knepley       suffix: hex_hash
1973285882aSMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_refine 1 -dm_plex_hash_location
198c4762a1bSJed Brown 
199*274abdb2SMatthew G. Knepley   testset:
200*274abdb2SMatthew G. Knepley     args: -centroids 0 -custom \
201*274abdb2SMatthew G. Knepley           -dm_plex_simplex 0 -dm_plex_box_faces 21,21 -dm_distribute_overlap 4 -petscpartitioner_type simple
202*274abdb2SMatthew G. Knepley     nsize: 2
203*274abdb2SMatthew G. Knepley 
204*274abdb2SMatthew G. Knepley     test:
205*274abdb2SMatthew G. Knepley       suffix: quad_overlap
206*274abdb2SMatthew G. Knepley       args: -dm_plex_hash_location {{0 1}}
207*274abdb2SMatthew G. Knepley 
208c4762a1bSJed Brown TEST*/
209