xref: /petsc/src/dm/impls/plex/tests/ex17.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   CHKERRQ(DMCreate(comm, dm));
10   CHKERRQ(DMSetType(*dm, DMPLEX));
11   CHKERRQ(DMSetFromOptions(*dm));
12   CHKERRQ(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   CHKERRQ(DMGetCoordinateDim(dm, &cdim));
27   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
28   /* Locate all centroids */
29   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart)*cdim, &points));
30   CHKERRQ(VecSetBlockSize(points, cdim));
31   CHKERRQ(VecGetArray(points, &a));
32   for (c = cStart; c < cEnd; ++c) {
33     PetscReal          centroid[3];
34     PetscInt           off = (c - cStart)*cdim, d;
35 
36     CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
37     for (d = 0; d < cdim; ++d) a[off+d] = centroid[d];
38   }
39   CHKERRQ(VecRestoreArray(points, &a));
40   CHKERRQ(DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF));
41   CHKERRQ(VecDestroy(&points));
42   CHKERRQ(PetscSFGetGraph(cellSF, NULL, &n, NULL, &cells));
43   if (n != (cEnd - cStart)) {
44     for (c = 0; c < n; ++c) {
45       if (cells[c].index != c+cStart) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Could not locate centroid of cell %D, error %D\n", c+cStart, cells[c].index));
46     }
47     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Located %D points instead of %D", n, cEnd - cStart);
48   }
49   for (c = cStart; c < cEnd; ++c) {
50     PetscCheckFalse(cells[c - cStart].index != c,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not locate centroid of cell %D, instead found %D", c, cells[c - cStart].index);
51   }
52   CHKERRQ(PetscSFDestroy(&cellSF));
53   PetscFunctionReturn(0);
54 }
55 
56 int main(int argc, char **argv)
57 {
58   DM             dm;
59 
60   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
61   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &dm));
62   CHKERRQ(TestLocation(dm));
63   CHKERRQ(DMDestroy(&dm));
64   CHKERRQ(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