xref: /petsc/src/dm/tutorials/ex20.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "DMSwarm-PIC demonstator of inserting points into a cell DM \n\
3c4762a1bSJed Brown Options: \n\
4c4762a1bSJed Brown -mode {0,1} : 0 ==> DMDA, 1 ==> DMPLEX cell DM \n\
5c4762a1bSJed Brown -dim {2,3}  : spatial dimension\n";
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petsc.h>
8c4762a1bSJed Brown #include <petscdm.h>
9c4762a1bSJed Brown #include <petscdmda.h>
10c4762a1bSJed Brown #include <petscdmplex.h>
11c4762a1bSJed Brown #include <petscdmswarm.h>
12c4762a1bSJed Brown 
139371c9d4SSatish Balay PetscErrorCode pic_insert_DMDA(PetscInt dim) {
14c4762a1bSJed Brown   DM        celldm = NULL, swarm;
15c4762a1bSJed Brown   PetscInt  dof, stencil_width;
16c4762a1bSJed Brown   PetscReal min[3], max[3];
17c4762a1bSJed Brown   PetscInt  ndir[3];
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   PetscFunctionBegin;
20c4762a1bSJed Brown   /* Create the background cell DM */
21c4762a1bSJed Brown   dof           = 1;
22c4762a1bSJed Brown   stencil_width = 1;
23*48a46eb9SPierre Jolivet   if (dim == 2) PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 25, 13, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, &celldm));
24*48a46eb9SPierre Jolivet   if (dim == 3) PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 25, 13, 19, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, stencil_width, NULL, NULL, NULL, &celldm));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(DMDASetElementType(celldm, DMDA_ELEMENT_Q1));
279566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(celldm));
289566063dSJacob Faibussowitsch   PetscCall(DMSetUp(celldm));
29c4762a1bSJed Brown 
309566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(celldm, 0.0, 2.0, 0.0, 1.0, 0.0, 1.5));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Create the DMSwarm */
339566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
349566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm"));
359566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm, DMSWARM));
369566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm, dim));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
399566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
409566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm, celldm));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
439566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE));
449566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE));
459566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
489566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
519566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_REGULAR, 3));
529371c9d4SSatish Balay   min[0]  = 0.5;
539371c9d4SSatish Balay   max[0]  = 0.7;
549371c9d4SSatish Balay   min[1]  = 0.5;
559371c9d4SSatish Balay   max[1]  = 0.8;
569371c9d4SSatish Balay   min[2]  = 0.5;
579371c9d4SSatish Balay   max[2]  = 0.9;
58c4762a1bSJed Brown   ndir[0] = ndir[1] = ndir[2] = 30;
599566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetPointsUniformCoordinates(swarm, min, max, ndir, ADD_VALUES));
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* This should be dispatched from a regular DMView() */
629566063dSJacob Faibussowitsch   PetscCall(DMSwarmViewXDMF(swarm, "ex20.xmf"));
639566063dSJacob Faibussowitsch   PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD));
649566063dSJacob Faibussowitsch   PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   {
67c4762a1bSJed Brown     PetscInt    npoints, *list;
68c4762a1bSJed Brown     PetscMPIInt rank;
69c4762a1bSJed Brown 
709566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
719566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetAccess(swarm));
729566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetNumberOfPointsPerCell(swarm, 0, &npoints));
739566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(swarm, rank, &npoints, &list));
749566063dSJacob Faibussowitsch     PetscCall(PetscFree(list));
759566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortRestoreAccess(swarm));
76c4762a1bSJed Brown   }
779566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate(swarm, PETSC_FALSE));
789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
799566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
80c4762a1bSJed Brown   PetscFunctionReturn(0);
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
839371c9d4SSatish Balay PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim) {
84c4762a1bSJed Brown   DM          celldm = NULL, swarm, distributedMesh = NULL;
85c4762a1bSJed Brown   const char *fieldnames[] = {"viscosity"};
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   PetscFunctionBegin;
88c4762a1bSJed Brown   /* Create the background cell DM */
89c4762a1bSJed Brown   if (dim == 2) {
90c4762a1bSJed Brown     PetscInt   cells_per_dim[2], nx[2];
91c4762a1bSJed Brown     PetscInt   n_tricells;
92c4762a1bSJed Brown     PetscInt   n_trivert;
93a4a685f2SJacob Faibussowitsch     PetscInt  *tricells;
94a4a685f2SJacob Faibussowitsch     PetscReal *trivert, dx, dy;
95c4762a1bSJed Brown     PetscInt   ii, jj, cnt;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown     cells_per_dim[0] = 4;
98c4762a1bSJed Brown     cells_per_dim[1] = 4;
99c4762a1bSJed Brown     n_tricells       = cells_per_dim[0] * cells_per_dim[1] * 2;
100c4762a1bSJed Brown     nx[0]            = cells_per_dim[0] + 1;
101c4762a1bSJed Brown     nx[1]            = cells_per_dim[1] + 1;
102c4762a1bSJed Brown     n_trivert        = nx[0] * nx[1];
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_tricells * 3, &tricells));
1059566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nx[0] * nx[1] * 2, &trivert));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown     /* verts */
108c4762a1bSJed Brown     cnt = 0;
109c4762a1bSJed Brown     dx  = 2.0 / ((PetscReal)cells_per_dim[0]);
110c4762a1bSJed Brown     dy  = 1.0 / ((PetscReal)cells_per_dim[1]);
111c4762a1bSJed Brown     for (jj = 0; jj < nx[1]; jj++) {
112c4762a1bSJed Brown       for (ii = 0; ii < nx[0]; ii++) {
113c4762a1bSJed Brown         trivert[2 * cnt + 0] = 0.0 + ii * dx;
114c4762a1bSJed Brown         trivert[2 * cnt + 1] = 0.0 + jj * dy;
115c4762a1bSJed Brown         cnt++;
116c4762a1bSJed Brown       }
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     /* connectivity */
120c4762a1bSJed Brown     cnt = 0;
121c4762a1bSJed Brown     for (jj = 0; jj < cells_per_dim[1]; jj++) {
122c4762a1bSJed Brown       for (ii = 0; ii < cells_per_dim[0]; ii++) {
123c4762a1bSJed Brown         PetscInt idx, idx0, idx1, idx2, idx3;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown         idx  = (ii) + (jj)*nx[0];
126c4762a1bSJed Brown         idx0 = idx;
127c4762a1bSJed Brown         idx1 = idx0 + 1;
128c4762a1bSJed Brown         idx2 = idx1 + nx[0];
129c4762a1bSJed Brown         idx3 = idx0 + nx[0];
130c4762a1bSJed Brown 
131c4762a1bSJed Brown         tricells[3 * cnt + 0] = idx0;
132c4762a1bSJed Brown         tricells[3 * cnt + 1] = idx1;
133c4762a1bSJed Brown         tricells[3 * cnt + 2] = idx2;
134c4762a1bSJed Brown         cnt++;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown         tricells[3 * cnt + 0] = idx0;
137c4762a1bSJed Brown         tricells[3 * cnt + 1] = idx2;
138c4762a1bSJed Brown         tricells[3 * cnt + 2] = idx3;
139c4762a1bSJed Brown         cnt++;
140c4762a1bSJed Brown       }
141c4762a1bSJed Brown     }
1429566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, n_tricells, n_trivert, 3, PETSC_TRUE, tricells, dim, trivert, &celldm));
1439566063dSJacob Faibussowitsch     PetscCall(PetscFree(trivert));
1449566063dSJacob Faibussowitsch     PetscCall(PetscFree(tricells));
145c4762a1bSJed Brown   }
14608401ef6SPierre Jolivet   PetscCheck(dim != 3, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Only 2D PLEX example supported");
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Distribute mesh over processes */
1499566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh));
150c4762a1bSJed Brown   if (distributedMesh) {
1519566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&celldm));
152c4762a1bSJed Brown     celldm = distributedMesh;
153c4762a1bSJed Brown   }
1549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)celldm, "Cells"));
1559566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(celldm));
156c4762a1bSJed Brown   {
157c4762a1bSJed Brown     PetscInt     numComp[] = {1};
158c4762a1bSJed Brown     PetscInt     numDof[]  = {1, 0, 0}; /* vert, edge, cell */
159c4762a1bSJed Brown     PetscInt     numBC     = 0;
160c4762a1bSJed Brown     PetscSection section;
161c4762a1bSJed Brown 
1629566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &section));
1639566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(celldm, section));
1649566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
165c4762a1bSJed Brown   }
1669566063dSJacob Faibussowitsch   PetscCall(DMSetUp(celldm));
167c4762a1bSJed Brown   {
168c4762a1bSJed Brown     PetscViewer viewer;
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
1719566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK));
1729566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
1739566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(viewer, "ex20plex.vtk"));
1749566063dSJacob Faibussowitsch     PetscCall(DMView(celldm, viewer));
1759566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
176c4762a1bSJed Brown   }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Create the DMSwarm */
1799566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm"));
1819566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm, DMSWARM));
1829566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm, dim));
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
1859566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm, celldm));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
1889566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE));
1899566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE));
1909566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
1939566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
1969566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_SUBDIVISION, 2));
1979566063dSJacob Faibussowitsch   PetscCall(DMSwarmViewFieldsXDMF(swarm, "ex20.xmf", 1, fieldnames));
1989566063dSJacob Faibussowitsch   PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD));
1999566063dSJacob Faibussowitsch   PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD));
2009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
2019566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
202c4762a1bSJed Brown   PetscFunctionReturn(0);
203c4762a1bSJed Brown }
204c4762a1bSJed Brown 
2059371c9d4SSatish Balay PetscErrorCode pic_insert_DMPLEX(PetscBool is_simplex, PetscInt dim) {
206c4762a1bSJed Brown   DM          celldm, swarm, distributedMesh = NULL;
207c4762a1bSJed Brown   const char *fieldnames[] = {"viscosity", "DMSwarm_rank"};
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBegin;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* Create the background cell DM */
212c4762a1bSJed Brown   {
213c4762a1bSJed Brown     PetscInt faces[3] = {4, 2, 4};
2149566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm));
215c4762a1bSJed Brown   }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* Distribute mesh over processes */
2189566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh));
219c4762a1bSJed Brown   if (distributedMesh) {
2209566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&celldm));
221c4762a1bSJed Brown     celldm = distributedMesh;
222c4762a1bSJed Brown   }
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)celldm, "Cells"));
2249566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(celldm));
225c4762a1bSJed Brown   {
226c4762a1bSJed Brown     PetscInt     numComp[] = {1};
227c4762a1bSJed Brown     PetscInt     numDof[]  = {1, 0, 0}; /* vert, edge, cell */
228c4762a1bSJed Brown     PetscInt     numBC     = 0;
229c4762a1bSJed Brown     PetscSection section;
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, &section));
2329566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(celldm, section));
2339566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
234c4762a1bSJed Brown   }
2359566063dSJacob Faibussowitsch   PetscCall(DMSetUp(celldm));
236c4762a1bSJed Brown   {
237c4762a1bSJed Brown     PetscViewer viewer;
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
2409566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK));
2419566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
2429566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(viewer, "ex20plex.vtk"));
2439566063dSJacob Faibussowitsch     PetscCall(DMView(celldm, viewer));
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
245c4762a1bSJed Brown   }
246c4762a1bSJed Brown 
2479566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
2489566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)swarm, "Swarm"));
2499566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm, DMSWARM));
2509566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm, dim));
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
2539566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm, celldm));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
2569566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "viscosity", 1, PETSC_DOUBLE));
2579566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "density", 1, PETSC_DOUBLE));
2589566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
2619566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
2649566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_GAUSS, 3));
2659566063dSJacob Faibussowitsch   PetscCall(DMSwarmViewFieldsXDMF(swarm, "ex20.xmf", 2, fieldnames));
2669566063dSJacob Faibussowitsch   PetscCall(DMView(celldm, PETSC_VIEWER_STDOUT_WORLD));
2679566063dSJacob Faibussowitsch   PetscCall(DMView(swarm, PETSC_VIEWER_STDOUT_WORLD));
2689566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
2699566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
270c4762a1bSJed Brown   PetscFunctionReturn(0);
271c4762a1bSJed Brown }
272c4762a1bSJed Brown 
2739371c9d4SSatish Balay int main(int argc, char **args) {
274c4762a1bSJed Brown   PetscInt mode = 0;
275c4762a1bSJed Brown   PetscInt dim  = 2;
276c4762a1bSJed Brown 
277327415f7SBarry Smith   PetscFunctionBeginUser;
2789566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
2799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mode", &mode, NULL));
2809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
281c4762a1bSJed Brown   switch (mode) {
2829371c9d4SSatish Balay   case 0: PetscCall(pic_insert_DMDA(dim)); break;
283c4762a1bSJed Brown   case 1:
284c4762a1bSJed Brown     /* tri / tet */
2859566063dSJacob Faibussowitsch     PetscCall(pic_insert_DMPLEX(PETSC_TRUE, dim));
286c4762a1bSJed Brown     break;
287c4762a1bSJed Brown   case 2:
288c4762a1bSJed Brown     /* quad / hex */
2899566063dSJacob Faibussowitsch     PetscCall(pic_insert_DMPLEX(PETSC_FALSE, dim));
290c4762a1bSJed Brown     break;
2919371c9d4SSatish Balay   default: PetscCall(pic_insert_DMDA(dim)); break;
292c4762a1bSJed Brown   }
2939566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
294b122ec5aSJacob Faibussowitsch   return 0;
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown /*TEST
298c4762a1bSJed Brown 
299c4762a1bSJed Brown    test:
300c4762a1bSJed Brown       args:
301c4762a1bSJed Brown       requires: !complex double
3027ed4f029SJed Brown       filter: grep -v atomic
303c4762a1bSJed Brown       filter_output: grep -v atomic
304c4762a1bSJed Brown 
305c4762a1bSJed Brown    test:
306c4762a1bSJed Brown       suffix: 2
307c4762a1bSJed Brown       requires: triangle double !complex
308c4762a1bSJed Brown       args: -mode 1
3097ed4f029SJed Brown       filter: grep -v atomic
310c4762a1bSJed Brown       filter_output: grep -v atomic
311c4762a1bSJed Brown 
312c4762a1bSJed Brown TEST*/
313