xref: /petsc/src/dm/impls/plex/tests/ex17.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests for point location\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscsf.h>
4*c4762a1bSJed Brown #include <petscdmplex.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown typedef struct {
7*c4762a1bSJed Brown   /* Domain and mesh definition */
8*c4762a1bSJed Brown   PetscInt  dim;                          /* The topological mesh dimension */
9*c4762a1bSJed Brown   PetscBool cellSimplex;                  /* Use simplices or hexes */
10*c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
11*c4762a1bSJed Brown   PetscBool testPartition;                /* Use a fixed partitioning for testing */
12*c4762a1bSJed Brown   PetscInt  testNum;                      /* Labels the different test partitions */
13*c4762a1bSJed Brown } AppCtx;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
16*c4762a1bSJed Brown {
17*c4762a1bSJed Brown   PetscErrorCode ierr;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   PetscFunctionBeginUser;
20*c4762a1bSJed Brown   options->dim           = 2;
21*c4762a1bSJed Brown   options->cellSimplex   = PETSC_TRUE;
22*c4762a1bSJed Brown   options->filename[0]   = '\0';
23*c4762a1bSJed Brown   options->testPartition = PETSC_TRUE;
24*c4762a1bSJed Brown   options->testNum       = 0;
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex13.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex13.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "The mesh file", "ex13.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex13.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-test_num", "The test partition number", "ex13.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
33*c4762a1bSJed Brown   PetscFunctionReturn(0);
34*c4762a1bSJed Brown }
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
37*c4762a1bSJed Brown {
38*c4762a1bSJed Brown   DM             dmDist      = NULL;
39*c4762a1bSJed Brown   PetscInt       dim         = user->dim;
40*c4762a1bSJed Brown   PetscBool      cellSimplex = user->cellSimplex;
41*c4762a1bSJed Brown   const char    *filename    = user->filename;
42*c4762a1bSJed Brown   size_t         len;
43*c4762a1bSJed Brown   PetscErrorCode ierr;
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown   PetscFunctionBeginUser;
46*c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
47*c4762a1bSJed Brown   if (len) {ierr = DMPlexCreateFromFile(comm, filename, PETSC_TRUE, dm);CHKERRQ(ierr);}
48*c4762a1bSJed Brown   else     {ierr = DMPlexCreateBoxMesh(comm, dim, cellSimplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);}
49*c4762a1bSJed Brown   if (user->testPartition) {
50*c4762a1bSJed Brown     PetscPartitioner part;
51*c4762a1bSJed Brown     PetscInt         *sizes  = NULL;
52*c4762a1bSJed Brown     PetscInt         *points = NULL;
53*c4762a1bSJed Brown     PetscMPIInt      rank, size;
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown     ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
56*c4762a1bSJed Brown     ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
57*c4762a1bSJed Brown     if (!rank) {
58*c4762a1bSJed Brown       if (dim == 2 && cellSimplex && size == 2) {
59*c4762a1bSJed Brown         switch (user->testNum) {
60*c4762a1bSJed Brown         case 0: {
61*c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {4, 4};
62*c4762a1bSJed Brown           PetscInt triPoints_p2[8] = {3, 5, 6, 7, 0, 1, 2, 4};
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr);
65*c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
66*c4762a1bSJed Brown           ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr);
67*c4762a1bSJed Brown           break;}
68*c4762a1bSJed Brown         case 1: {
69*c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {6, 2};
70*c4762a1bSJed Brown           PetscInt triPoints_p2[8] = {1, 2, 3, 4, 6, 7, 0, 5};
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 8, &points);CHKERRQ(ierr);
73*c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
74*c4762a1bSJed Brown           ierr = PetscArraycpy(points, triPoints_p2, 8);CHKERRQ(ierr);
75*c4762a1bSJed Brown           break;}
76*c4762a1bSJed Brown         default:
77*c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum);
78*c4762a1bSJed Brown         }
79*c4762a1bSJed Brown       } else if (dim == 2 && cellSimplex && size == 3) {
80*c4762a1bSJed Brown         PetscInt triSizes_p3[3]  = {3, 3, 2};
81*c4762a1bSJed Brown         PetscInt triPoints_p3[8] = {1, 2, 4, 3, 6, 7, 0, 5};
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown         ierr = PetscMalloc2(3, &sizes, 8, &points);CHKERRQ(ierr);
84*c4762a1bSJed Brown         ierr = PetscArraycpy(sizes,  triSizes_p3, 3);CHKERRQ(ierr);
85*c4762a1bSJed Brown         ierr = PetscArraycpy(points, triPoints_p3, 8);CHKERRQ(ierr);
86*c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && size == 2) {
87*c4762a1bSJed Brown         PetscInt quadSizes_p2[2]  = {2, 2};
88*c4762a1bSJed Brown         PetscInt quadPoints_p2[4] = {2, 3, 0, 1};
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown         ierr = PetscMalloc2(2, &sizes, 4, &points);CHKERRQ(ierr);
91*c4762a1bSJed Brown         ierr = PetscArraycpy(sizes,  quadSizes_p2, 2);CHKERRQ(ierr);
92*c4762a1bSJed Brown         ierr = PetscArraycpy(points, quadPoints_p2, 4);CHKERRQ(ierr);
93*c4762a1bSJed Brown       } else SETERRQ3(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition dim: %D simplex: %D size: %D", dim, (PetscInt) cellSimplex, (PetscInt) size);
94*c4762a1bSJed Brown     }
95*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
96*c4762a1bSJed Brown     ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr);
97*c4762a1bSJed Brown     ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr);
98*c4762a1bSJed Brown     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
99*c4762a1bSJed Brown   }
100*c4762a1bSJed Brown   ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr);
101*c4762a1bSJed Brown   if (dmDist) {
102*c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
103*c4762a1bSJed Brown     *dm  = dmDist;
104*c4762a1bSJed Brown   }
105*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, cellSimplex ? "Simplicial Mesh" : "Tensor Product Mesh");CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
107*c4762a1bSJed Brown   PetscFunctionReturn(0);
108*c4762a1bSJed Brown }
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown static PetscErrorCode TestLocation(DM dm, AppCtx *user)
111*c4762a1bSJed Brown {
112*c4762a1bSJed Brown   PetscInt       dim = user->dim;
113*c4762a1bSJed Brown   PetscInt       cStart, cEnd, c;
114*c4762a1bSJed Brown   PetscErrorCode ierr;
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   PetscFunctionBeginUser;
117*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
118*c4762a1bSJed Brown   /* Locate all centroids */
119*c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
120*c4762a1bSJed Brown     Vec                v;
121*c4762a1bSJed Brown     PetscSF            cellSF = NULL;
122*c4762a1bSJed Brown     const PetscSFNode *cells;
123*c4762a1bSJed Brown     PetscScalar       *a;
124*c4762a1bSJed Brown     PetscReal          centroid[3];
125*c4762a1bSJed Brown     PetscInt           n, d;
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown     ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
128*c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF, dim, &v);CHKERRQ(ierr);
129*c4762a1bSJed Brown     ierr = VecSetBlockSize(v, dim);CHKERRQ(ierr);
130*c4762a1bSJed Brown     ierr = VecGetArray(v, &a);CHKERRQ(ierr);
131*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) a[d] = centroid[d];
132*c4762a1bSJed Brown     ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
133*c4762a1bSJed Brown     ierr = DMLocatePoints(dm, v, DM_POINTLOCATION_NONE, &cellSF);CHKERRQ(ierr);
134*c4762a1bSJed Brown     ierr = VecDestroy(&v);CHKERRQ(ierr);
135*c4762a1bSJed Brown     ierr = PetscSFGetGraph(cellSF,NULL,&n,NULL,&cells);CHKERRQ(ierr);
136*c4762a1bSJed Brown     if (n        != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Found %d cells instead %d", n, 1);
137*c4762a1bSJed Brown     if (cells[0].index != c) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not locate centroid of cell %d, instead found %d", c, cells[0].index);
138*c4762a1bSJed Brown     ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
139*c4762a1bSJed Brown   }
140*c4762a1bSJed Brown   PetscFunctionReturn(0);
141*c4762a1bSJed Brown }
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown int main(int argc, char **argv)
144*c4762a1bSJed Brown {
145*c4762a1bSJed Brown   DM             dm;
146*c4762a1bSJed Brown   AppCtx         user; /* user-defined work context */
147*c4762a1bSJed Brown   PetscErrorCode ierr;
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
150*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = TestLocation(dm, &user);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = PetscFinalize();
155*c4762a1bSJed Brown   return ierr;
156*c4762a1bSJed Brown }
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown /*TEST
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown   test:
161*c4762a1bSJed Brown     suffix: 0
162*c4762a1bSJed Brown     requires: triangle
163*c4762a1bSJed Brown     args: -test_partition 0 -dm_view ascii::ascii_info_detail
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown TEST*/
166