xref: /petsc/src/dm/impls/plex/tests/ex17.c (revision 3285882ab300fffbe2ae6c3a9793ee7447c414ab)
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 
630602db0SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   PetscErrorCode ierr;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   PetscFunctionBeginUser;
1130602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1230602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1330602db0SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
14c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
15c4762a1bSJed Brown   PetscFunctionReturn(0);
16c4762a1bSJed Brown }
17c4762a1bSJed Brown 
1830602db0SMatthew G. Knepley static PetscErrorCode TestLocation(DM dm)
19c4762a1bSJed Brown {
20*3285882aSMatthew G. Knepley   Vec                points;
21*3285882aSMatthew G. Knepley   PetscSF            cellSF = NULL;
22*3285882aSMatthew G. Knepley   const PetscSFNode *cells;
23*3285882aSMatthew G. Knepley   PetscScalar       *a;
24*3285882aSMatthew G. Knepley   PetscInt           cdim, n;
25c4762a1bSJed Brown   PetscInt           cStart, cEnd, c;
26c4762a1bSJed Brown   PetscErrorCode     ierr;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   PetscFunctionBeginUser;
29*3285882aSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
31c4762a1bSJed Brown   /* Locate all centroids */
32*3285882aSMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, (cEnd - cStart)*cdim, &points);CHKERRQ(ierr);
33*3285882aSMatthew G. Knepley   ierr = VecSetBlockSize(points, cdim);CHKERRQ(ierr);
34*3285882aSMatthew G. Knepley   ierr = VecGetArray(points, &a);CHKERRQ(ierr);
35c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
36c4762a1bSJed Brown     PetscReal          centroid[3];
37*3285882aSMatthew G. Knepley     PetscInt           off = (c - cStart)*cdim, d;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown     ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
40*3285882aSMatthew G. Knepley     for (d = 0; d < cdim; ++d) a[off+d] = centroid[d];
41c4762a1bSJed Brown   }
42*3285882aSMatthew G. Knepley   ierr = VecRestoreArray(points, &a);CHKERRQ(ierr);
43*3285882aSMatthew G. Knepley   ierr = DMLocatePoints(dm, points, DM_POINTLOCATION_NONE, &cellSF);CHKERRQ(ierr);
44*3285882aSMatthew G. Knepley   ierr = VecDestroy(&points);CHKERRQ(ierr);
45*3285882aSMatthew G. Knepley   ierr = PetscSFGetGraph(cellSF, NULL, &n, NULL, &cells);CHKERRQ(ierr);
46*3285882aSMatthew G. Knepley   if (n != (cEnd - cStart)) {
47*3285882aSMatthew G. Knepley     for (c = 0; c < n; ++c) {
48*3285882aSMatthew G. Knepley       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*3285882aSMatthew G. Knepley     }
50*3285882aSMatthew G. Knepley     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Located %D points instead of %D", n, cEnd - cStart);
51*3285882aSMatthew G. Knepley   }
52*3285882aSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
53*3285882aSMatthew G. Knepley     if (cells[c - cStart].index != c) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not locate centroid of cell %D, instead found %D", c, cells[c - cStart].index);
54*3285882aSMatthew G. Knepley   }
55*3285882aSMatthew G. Knepley   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown int main(int argc, char **argv)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   DM             dm;
62c4762a1bSJed Brown   PetscErrorCode ierr;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
6530602db0SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);
6630602db0SMatthew G. Knepley   ierr = TestLocation(dm);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscFinalize();
69c4762a1bSJed Brown   return ierr;
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown /*TEST
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   test:
751af33867SMatthew G. Knepley     suffix: seg
761af33867SMatthew G. Knepley     args: -dm_plex_dim 1 -dm_plex_box_faces 10
771af33867SMatthew G. Knepley 
78*3285882aSMatthew G. Knepley   testset:
79*3285882aSMatthew G. Knepley     args: -dm_plex_box_faces 5,5
80*3285882aSMatthew G. Knepley 
811af33867SMatthew G. Knepley     test:
821af33867SMatthew G. Knepley       suffix: tri
83c4762a1bSJed Brown       requires: triangle
84*3285882aSMatthew G. Knepley 
85*3285882aSMatthew G. Knepley     test:
86*3285882aSMatthew G. Knepley       suffix: tri_hash
87*3285882aSMatthew G. Knepley       requires: triangle
88*3285882aSMatthew G. Knepley       args: -dm_refine 2 -dm_plex_hash_location
891af33867SMatthew G. Knepley 
901af33867SMatthew G. Knepley     test:
911af33867SMatthew G. Knepley       suffix: quad
92*3285882aSMatthew G. Knepley       args: -dm_plex_simplex 0
93*3285882aSMatthew G. Knepley 
94*3285882aSMatthew G. Knepley     test:
95*3285882aSMatthew G. Knepley       suffix: quad_hash
96*3285882aSMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_refine 2 -dm_plex_hash_location
97*3285882aSMatthew G. Knepley 
98*3285882aSMatthew G. Knepley   testset:
99*3285882aSMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3
1001af33867SMatthew G. Knepley 
1011af33867SMatthew G. Knepley     test:
1021af33867SMatthew G. Knepley       suffix: tet
1031af33867SMatthew G. Knepley       requires: ctetgen
104*3285882aSMatthew G. Knepley 
105*3285882aSMatthew G. Knepley     test:
106*3285882aSMatthew G. Knepley       suffix: tet_hash
107*3285882aSMatthew G. Knepley       requires: ctetgen
108*3285882aSMatthew G. Knepley       args: -dm_refine 1 -dm_plex_hash_location
1091af33867SMatthew G. Knepley 
1101af33867SMatthew G. Knepley     test:
1111af33867SMatthew G. Knepley       suffix: hex
112*3285882aSMatthew G. Knepley       args: -dm_plex_simplex 0
113*3285882aSMatthew G. Knepley 
114*3285882aSMatthew G. Knepley     test:
115*3285882aSMatthew G. Knepley       suffix: hex_hash
116*3285882aSMatthew G. Knepley       args: -dm_plex_simplex 0 -dm_refine 1 -dm_plex_hash_location
117c4762a1bSJed Brown 
118c4762a1bSJed Brown TEST*/
119