xref: /petsc/src/dm/impls/plex/tests/ex17.c (revision 217f7fc6b7ccbb041bf230fae922d6b0eefae257)
1 static char help[] = "Tests for point location\n\n";
2 
3 #include <petscsf.h>
4 #include <petscdmplex.h>
5 
6 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
7 {
8   PetscFunctionBeginUser;
9   PetscCall(DMCreate(comm, dm));
10   PetscCall(DMSetType(*dm, DMPLEX));
11   PetscCall(DMSetFromOptions(*dm));
12   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
13   PetscFunctionReturn(0);
14 }
15 
16 static PetscErrorCode TestLocation(DM dm)
17 {
18   Vec                points;
19   PetscSF            cellSF = NULL;
20   const PetscSFNode *cells;
21   PetscScalar       *a;
22   PetscInt           cdim, n;
23   PetscInt           cStart, cEnd, c;
24 
25   PetscFunctionBeginUser;
26   PetscCall(DMGetCoordinateDim(dm, &cdim));
27   PetscCall(DMGetCoordinatesLocalSetUp(dm));
28   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
29   /* Locate all centroids */
30   PetscCall(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart) * cdim, &points));
31   PetscCall(VecSetBlockSize(points, cdim));
32   PetscCall(VecGetArray(points, &a));
33   for (c = cStart; c < cEnd; ++c) {
34     PetscReal centroid[3];
35     PetscInt  off = (c - cStart) * cdim, d;
36 
37     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
38     for (d = 0; d < cdim; ++d) a[off + d] = centroid[d];
39   }
40   PetscCall(VecRestoreArray(points, &a));
41   PetscCall(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF));
42   PetscCall(VecDestroy(&points));
43   PetscCall(PetscSFGetGraph(cellSF, NULL, &n, NULL, &cells));
44   if (n != (cEnd - cStart)) {
45     for (c = 0; c < n; ++c) {
46       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));
47     }
48     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Located %" PetscInt_FMT " points instead of %" PetscInt_FMT, n, cEnd - cStart);
49   }
50   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);
51   PetscCall(PetscSFDestroy(&cellSF));
52   PetscFunctionReturn(0);
53 }
54 
55 int main(int argc, char **argv)
56 {
57   DM dm;
58 
59   PetscFunctionBeginUser;
60   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
61   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
62   PetscCall(TestLocation(dm));
63   PetscCall(DMDestroy(&dm));
64   PetscCall(PetscFinalize());
65   return 0;
66 }
67 
68 /*TEST
69 
70   testset:
71     args: -dm_plex_dim 1 -dm_plex_box_faces 10
72 
73     test:
74       suffix: seg
75 
76     test:
77       suffix: seg_hash
78       args: -dm_refine 2 -dm_plex_hash_location
79 
80   testset:
81     args: -dm_plex_box_faces 5,5
82 
83     test:
84       suffix: tri
85       requires: triangle
86 
87     test:
88       suffix: tri_hash
89       requires: triangle
90       args: -dm_refine 2 -dm_plex_hash_location
91 
92     test:
93       suffix: quad
94       args: -dm_plex_simplex 0
95 
96     test:
97       suffix: quad_hash
98       args: -dm_plex_simplex 0 -dm_refine 2 -dm_plex_hash_location
99 
100   testset:
101     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3
102 
103     test:
104       suffix: tet
105       requires: ctetgen
106 
107     test:
108       suffix: tet_hash
109       requires: ctetgen
110       args: -dm_refine 1 -dm_plex_hash_location
111 
112     test:
113       suffix: hex
114       args: -dm_plex_simplex 0
115 
116     test:
117       suffix: hex_hash
118       args: -dm_plex_simplex 0 -dm_refine 1 -dm_plex_hash_location
119 
120 TEST*/
121