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