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
665f07beaSMatthew G. Knepley /* To inspect the location process, use
765f07beaSMatthew G. Knepley
865f07beaSMatthew G. Knepley -dm_plex_print_locate 5 -info :dm
965f07beaSMatthew G. Knepley */
1065f07beaSMatthew G. Knepley
11274abdb2SMatthew G. Knepley typedef struct {
12274abdb2SMatthew G. Knepley PetscBool centroids;
13274abdb2SMatthew G. Knepley PetscBool custom;
14274abdb2SMatthew G. Knepley } AppCtx;
15274abdb2SMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)16274abdb2SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
17274abdb2SMatthew G. Knepley {
18274abdb2SMatthew G. Knepley PetscFunctionBeginUser;
19274abdb2SMatthew G. Knepley options->centroids = PETSC_TRUE;
20274abdb2SMatthew G. Knepley options->custom = PETSC_FALSE;
21274abdb2SMatthew G. Knepley
22274abdb2SMatthew G. Knepley PetscOptionsBegin(comm, "", "Point Location Options", "DMPLEX");
23274abdb2SMatthew G. Knepley PetscCall(PetscOptionsBool("-centroids", "Locate cell centroids", "ex17.c", options->centroids, &options->centroids, NULL));
24274abdb2SMatthew G. Knepley PetscCall(PetscOptionsBool("-custom", "Locate user-defined points", "ex17.c", options->custom, &options->custom, NULL));
25274abdb2SMatthew G. Knepley PetscOptionsEnd();
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27274abdb2SMatthew G. Knepley }
28274abdb2SMatthew G. Knepley
CreateMesh(MPI_Comm comm,DM * dm)29d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown PetscFunctionBeginUser;
329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37c4762a1bSJed Brown }
38c4762a1bSJed Brown
TestCentroidLocation(DM dm,AppCtx * user)39274abdb2SMatthew G. Knepley static PetscErrorCode TestCentroidLocation(DM dm, AppCtx *user)
40d71ae5a4SJacob Faibussowitsch {
413285882aSMatthew G. Knepley Vec points;
423285882aSMatthew G. Knepley PetscSF cellSF = NULL;
433285882aSMatthew G. Knepley const PetscSFNode *cells;
443285882aSMatthew G. Knepley PetscScalar *a;
453285882aSMatthew G. Knepley PetscInt cdim, n;
46c4762a1bSJed Brown PetscInt cStart, cEnd, c;
47c4762a1bSJed Brown
48c4762a1bSJed Brown PetscFunctionBeginUser;
493ba16761SJacob Faibussowitsch if (!user->centroids) PetscFunctionReturn(PETSC_SUCCESS);
509566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim));
518fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm));
529566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
53c4762a1bSJed Brown /* Locate all centroids */
549566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart) * cdim, &points));
559566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(points, cdim));
569566063dSJacob Faibussowitsch PetscCall(VecGetArray(points, &a));
57c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) {
58c4762a1bSJed Brown PetscReal centroid[3];
593285882aSMatthew G. Knepley PetscInt off = (c - cStart) * cdim, d;
60c4762a1bSJed Brown
619566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
623285882aSMatthew G. Knepley for (d = 0; d < cdim; ++d) a[off + d] = centroid[d];
63c4762a1bSJed Brown }
649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(points, &a));
659566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF));
669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&points));
679566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(cellSF, NULL, &n, NULL, &cells));
683285882aSMatthew G. Knepley if (n != (cEnd - cStart)) {
693285882aSMatthew G. Knepley for (c = 0; c < n; ++c) {
7063a3b9bcSJacob 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));
713285882aSMatthew G. Knepley }
7263a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Located %" PetscInt_FMT " points instead of %" PetscInt_FMT, n, cEnd - cStart);
733285882aSMatthew G. Knepley }
74ad540459SPierre 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);
759566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&cellSF));
763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown
TestCustomLocation(DM dm,AppCtx * user)79274abdb2SMatthew G. Knepley static PetscErrorCode TestCustomLocation(DM dm, AppCtx *user)
80274abdb2SMatthew G. Knepley {
81274abdb2SMatthew G. Knepley PetscSF cellSF = NULL;
82274abdb2SMatthew G. Knepley const PetscSFNode *cells;
83274abdb2SMatthew G. Knepley const PetscInt *found;
84274abdb2SMatthew G. Knepley Vec points;
85274abdb2SMatthew G. Knepley PetscScalar coords[2] = {0.5, 0.5};
86274abdb2SMatthew G. Knepley PetscInt cdim, Np = 1, Nfd;
87274abdb2SMatthew G. Knepley PetscMPIInt rank;
88274abdb2SMatthew G. Knepley MPI_Comm comm;
89274abdb2SMatthew G. Knepley
90274abdb2SMatthew G. Knepley PetscFunctionBeginUser;
913ba16761SJacob Faibussowitsch if (!user->custom) PetscFunctionReturn(PETSC_SUCCESS);
92274abdb2SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim));
93274abdb2SMatthew G. Knepley
94274abdb2SMatthew G. Knepley // Locate serially on each process
95274abdb2SMatthew G. Knepley PetscCall(VecCreate(PETSC_COMM_SELF, &points));
96274abdb2SMatthew G. Knepley PetscCall(VecSetBlockSize(points, cdim));
97274abdb2SMatthew G. Knepley PetscCall(VecSetSizes(points, Np * cdim, PETSC_DETERMINE));
98274abdb2SMatthew G. Knepley PetscCall(VecSetFromOptions(points));
99274abdb2SMatthew G. Knepley for (PetscInt p = 0; p < Np; ++p) {
100274abdb2SMatthew G. Knepley const PetscInt idx[2] = {p * cdim, p * cdim + 1};
101274abdb2SMatthew G. Knepley PetscCall(VecSetValues(points, cdim, idx, coords, INSERT_VALUES));
102274abdb2SMatthew G. Knepley }
103274abdb2SMatthew G. Knepley PetscCall(VecAssemblyBegin(points));
104274abdb2SMatthew G. Knepley PetscCall(VecAssemblyEnd(points));
105274abdb2SMatthew G. Knepley
106274abdb2SMatthew G. Knepley PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF));
107274abdb2SMatthew G. Knepley
108274abdb2SMatthew G. Knepley PetscCall(PetscSFGetGraph(cellSF, NULL, &Nfd, &found, &cells));
109274abdb2SMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
110274abdb2SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(comm, &rank));
111274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, "[%d] Found %" PetscInt_FMT " particles\n", rank, Nfd));
112274abdb2SMatthew G. Knepley for (PetscInt p = 0; p < Nfd; ++p) {
113274abdb2SMatthew G. Knepley const PetscInt point = found ? found[p] : p;
114274abdb2SMatthew G. Knepley const PetscScalar *array;
115274abdb2SMatthew G. Knepley PetscScalar *ccoords = NULL;
116274abdb2SMatthew G. Knepley PetscInt numCoords;
117274abdb2SMatthew G. Knepley PetscBool isDG;
118274abdb2SMatthew G. Knepley
119274abdb2SMatthew G. Knepley // Since the v comm is SELF, rank is always 0
120274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, " point %" PetscInt_FMT " cell %" PetscInt_FMT "\n", point, cells[p].index));
121274abdb2SMatthew G. Knepley PetscCall(DMPlexGetCellCoordinates(dm, cells[p].index, &isDG, &numCoords, &array, &ccoords));
122274abdb2SMatthew G. Knepley for (PetscInt c = 0; c < numCoords / cdim; ++c) {
123274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, " "));
124274abdb2SMatthew G. Knepley for (PetscInt d = 0; d < cdim; ++d) PetscCall(PetscSynchronizedPrintf(comm, " %g", (double)PetscRealPart(ccoords[c * cdim + d])));
125274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedPrintf(comm, "\n"));
126274abdb2SMatthew G. Knepley }
127274abdb2SMatthew G. Knepley PetscCall(DMPlexRestoreCellCoordinates(dm, cells[p].index, &isDG, &numCoords, &array, &ccoords));
128274abdb2SMatthew G. Knepley }
129274abdb2SMatthew G. Knepley PetscCall(PetscSynchronizedFlush(comm, PETSC_STDOUT));
130274abdb2SMatthew G. Knepley
131274abdb2SMatthew G. Knepley PetscCall(PetscSFDestroy(&cellSF));
132274abdb2SMatthew G. Knepley PetscCall(VecDestroy(&points));
1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
134274abdb2SMatthew G. Knepley }
135274abdb2SMatthew G. Knepley
main(int argc,char ** argv)136d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
137d71ae5a4SJacob Faibussowitsch {
138c4762a1bSJed Brown DM dm;
139274abdb2SMatthew G. Knepley AppCtx user;
140c4762a1bSJed Brown
141327415f7SBarry Smith PetscFunctionBeginUser;
1429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
143274abdb2SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1449566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
145274abdb2SMatthew G. Knepley PetscCall(TestCentroidLocation(dm, &user));
146274abdb2SMatthew G. Knepley PetscCall(TestCustomLocation(dm, &user));
1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
1489566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
149b122ec5aSJacob Faibussowitsch return 0;
150c4762a1bSJed Brown }
151c4762a1bSJed Brown
152c4762a1bSJed Brown /*TEST
153c4762a1bSJed Brown
154b26b5bf9SMatthew G. Knepley testset:
155b26b5bf9SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 10
156*3886731fSPierre Jolivet output_file: output/empty.out
157b26b5bf9SMatthew G. Knepley
158c4762a1bSJed Brown test:
1591af33867SMatthew G. Knepley suffix: seg
160b26b5bf9SMatthew G. Knepley
161b26b5bf9SMatthew G. Knepley test:
162b26b5bf9SMatthew G. Knepley suffix: seg_hash
163b26b5bf9SMatthew G. Knepley args: -dm_refine 2 -dm_plex_hash_location
1641af33867SMatthew G. Knepley
1653285882aSMatthew G. Knepley testset:
1663285882aSMatthew G. Knepley args: -dm_plex_box_faces 5,5
167*3886731fSPierre Jolivet output_file: output/empty.out
1683285882aSMatthew G. Knepley
1691af33867SMatthew G. Knepley test:
1701af33867SMatthew G. Knepley suffix: tri
171c4762a1bSJed Brown requires: triangle
1723285882aSMatthew G. Knepley
1733285882aSMatthew G. Knepley test:
1743285882aSMatthew G. Knepley suffix: tri_hash
1753285882aSMatthew G. Knepley requires: triangle
1763285882aSMatthew G. Knepley args: -dm_refine 2 -dm_plex_hash_location
1771af33867SMatthew G. Knepley
1781af33867SMatthew G. Knepley test:
1791af33867SMatthew G. Knepley suffix: quad
1803285882aSMatthew G. Knepley args: -dm_plex_simplex 0
1813285882aSMatthew G. Knepley
1823285882aSMatthew G. Knepley test:
183dd301514SZach Atkins suffix: quad_order_2
184dd301514SZach Atkins args: -dm_plex_simplex 0 -dm_coord_petscspace_degree 2
185dd301514SZach Atkins
186dd301514SZach Atkins test:
1873285882aSMatthew G. Knepley suffix: quad_hash
1883285882aSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 2 -dm_plex_hash_location
1893285882aSMatthew G. Knepley
1903285882aSMatthew G. Knepley testset:
1913285882aSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3
192*3886731fSPierre Jolivet output_file: output/empty.out
1931af33867SMatthew G. Knepley
1941af33867SMatthew G. Knepley test:
1951af33867SMatthew G. Knepley suffix: tet
1961af33867SMatthew G. Knepley requires: ctetgen
1973285882aSMatthew G. Knepley
1983285882aSMatthew G. Knepley test:
1993285882aSMatthew G. Knepley suffix: tet_hash
2003285882aSMatthew G. Knepley requires: ctetgen
2013285882aSMatthew G. Knepley args: -dm_refine 1 -dm_plex_hash_location
2021af33867SMatthew G. Knepley
2031af33867SMatthew G. Knepley test:
2041af33867SMatthew G. Knepley suffix: hex
2053285882aSMatthew G. Knepley args: -dm_plex_simplex 0
2063285882aSMatthew G. Knepley
2073285882aSMatthew G. Knepley test:
2083285882aSMatthew G. Knepley suffix: hex_hash
2093285882aSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 1 -dm_plex_hash_location
210c4762a1bSJed Brown
211dd301514SZach Atkins test:
212dd301514SZach Atkins suffix: hex_order_2
213dd301514SZach Atkins args: -dm_plex_simplex 0 -dm_refine 1 -dm_coord_petscspace_degree 2
214dd301514SZach Atkins nsize: 2
215dd301514SZach Atkins
216dd301514SZach Atkins test:
217dd301514SZach Atkins suffix: hex_order_3
218dd301514SZach Atkins args: -dm_plex_simplex 0 -dm_coord_petscspace_degree 3
219dd301514SZach Atkins
220274abdb2SMatthew G. Knepley testset:
221274abdb2SMatthew G. Knepley args: -centroids 0 -custom \
222274abdb2SMatthew G. Knepley -dm_plex_simplex 0 -dm_plex_box_faces 21,21 -dm_distribute_overlap 4 -petscpartitioner_type simple
223274abdb2SMatthew G. Knepley nsize: 2
224274abdb2SMatthew G. Knepley
225274abdb2SMatthew G. Knepley test:
226274abdb2SMatthew G. Knepley suffix: quad_overlap
227274abdb2SMatthew G. Knepley args: -dm_plex_hash_location {{0 1}}
228274abdb2SMatthew G. Knepley
22965f07beaSMatthew G. Knepley # Test location on a Monge Manifold
23065f07beaSMatthew G. Knepley testset:
23165f07beaSMatthew G. Knepley args: -dm_refine 3 -dm_coord_space 0 \
23265f07beaSMatthew G. Knepley -dm_plex_option_phases proj_ -cdm_proj_dm_plex_coordinate_dim 3 -proj_dm_coord_space \
23365f07beaSMatthew G. Knepley -proj_dm_coord_remap -proj_dm_coord_map sinusoid -proj_dm_coord_map_params 0.1,1.,1.
234*3886731fSPierre Jolivet output_file: output/empty.out
23565f07beaSMatthew G. Knepley
23665f07beaSMatthew G. Knepley test:
23765f07beaSMatthew G. Knepley requires: triangle
23865f07beaSMatthew G. Knepley suffix: tri_monge
23965f07beaSMatthew G. Knepley
24065f07beaSMatthew G. Knepley test:
24165f07beaSMatthew G. Knepley suffix: quad_monge
24265f07beaSMatthew G. Knepley args: -dm_plex_simplex 0
24365f07beaSMatthew G. Knepley
244c4762a1bSJed Brown TEST*/
245