xref: /petsc/src/dm/tutorials/ex20.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 
13c4762a1bSJed Brown PetscErrorCode pic_insert_DMDA(PetscInt dim)
14c4762a1bSJed Brown {
15c4762a1bSJed Brown   DM             celldm = NULL,swarm;
16c4762a1bSJed Brown   PetscInt       dof,stencil_width;
17c4762a1bSJed Brown   PetscReal      min[3],max[3];
18c4762a1bSJed Brown   PetscInt       ndir[3];
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBegin;
21c4762a1bSJed Brown   /* Create the background cell DM */
22c4762a1bSJed Brown   dof = 1;
23c4762a1bSJed Brown   stencil_width = 1;
24c4762a1bSJed Brown   if (dim == 2) {
259566063dSJacob Faibussowitsch     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));
26c4762a1bSJed Brown   }
27c4762a1bSJed Brown   if (dim == 3) {
289566063dSJacob Faibussowitsch     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));
29c4762a1bSJed Brown   }
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(DMDASetElementType(celldm,DMDA_ELEMENT_Q1));
329566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(celldm));
339566063dSJacob Faibussowitsch   PetscCall(DMSetUp(celldm));
34c4762a1bSJed Brown 
359566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(celldm,0.0,2.0,0.0,1.0,0.0,1.5));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Create the DMSwarm */
389566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm));
399566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)swarm,"Swarm"));
409566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm,DMSWARM));
419566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm,dim));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
449566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC));
459566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm,celldm));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
489566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE));
499566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE));
509566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
539566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm,4,0));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
569566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,3));
57c4762a1bSJed Brown   min[0] = 0.5; max[0] = 0.7;
58c4762a1bSJed Brown   min[1] = 0.5; max[1] = 0.8;
59c4762a1bSJed Brown   min[2] = 0.5; max[2] = 0.9;
60c4762a1bSJed Brown   ndir[0] = ndir[1] = ndir[2] = 30;
619566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetPointsUniformCoordinates(swarm,min,max,ndir,ADD_VALUES));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* This should be dispatched from a regular DMView() */
649566063dSJacob Faibussowitsch   PetscCall(DMSwarmViewXDMF(swarm,"ex20.xmf"));
659566063dSJacob Faibussowitsch   PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
669566063dSJacob Faibussowitsch   PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   {
69c4762a1bSJed Brown     PetscInt    npoints,*list;
70c4762a1bSJed Brown     PetscMPIInt rank;
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
739566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetAccess(swarm));
749566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetNumberOfPointsPerCell(swarm,0,&npoints));
759566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(swarm,rank,&npoints,&list));
769566063dSJacob Faibussowitsch     PetscCall(PetscFree(list));
779566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortRestoreAccess(swarm));
78c4762a1bSJed Brown   }
799566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate(swarm,PETSC_FALSE));
809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
82c4762a1bSJed Brown   PetscFunctionReturn(0);
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim)
86c4762a1bSJed Brown {
87c4762a1bSJed Brown   DM             celldm = NULL,swarm,distributedMesh = NULL;
88c4762a1bSJed Brown   const  char    *fieldnames[] = {"viscosity"};
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   PetscFunctionBegin;
91c4762a1bSJed Brown   /* Create the background cell DM */
92c4762a1bSJed Brown   if (dim == 2) {
93c4762a1bSJed Brown     PetscInt   cells_per_dim[2],nx[2];
94c4762a1bSJed Brown     PetscInt   n_tricells;
95c4762a1bSJed Brown     PetscInt   n_trivert;
96a4a685f2SJacob Faibussowitsch     PetscInt   *tricells;
97a4a685f2SJacob Faibussowitsch     PetscReal  *trivert,dx,dy;
98c4762a1bSJed Brown     PetscInt   ii,jj,cnt;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown     cells_per_dim[0] = 4;
101c4762a1bSJed Brown     cells_per_dim[1] = 4;
102c4762a1bSJed Brown     n_tricells = cells_per_dim[0] * cells_per_dim[1] * 2;
103c4762a1bSJed Brown     nx[0] = cells_per_dim[0] + 1;
104c4762a1bSJed Brown     nx[1] = cells_per_dim[1] + 1;
105c4762a1bSJed Brown     n_trivert = nx[0] * nx[1];
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n_tricells*3,&tricells));
1089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nx[0]*nx[1]*2,&trivert));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown     /* verts */
111c4762a1bSJed Brown     cnt = 0;
112c4762a1bSJed Brown     dx = 2.0/((PetscReal)cells_per_dim[0]);
113c4762a1bSJed Brown     dy = 1.0/((PetscReal)cells_per_dim[1]);
114c4762a1bSJed Brown     for (jj=0; jj<nx[1]; jj++) {
115c4762a1bSJed Brown       for (ii=0; ii<nx[0]; ii++) {
116c4762a1bSJed Brown         trivert[2*cnt+0] = 0.0 + ii * dx;
117c4762a1bSJed Brown         trivert[2*cnt+1] = 0.0 + jj * dy;
118c4762a1bSJed Brown         cnt++;
119c4762a1bSJed Brown       }
120c4762a1bSJed Brown     }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown     /* connectivity */
123c4762a1bSJed Brown     cnt = 0;
124c4762a1bSJed Brown     for (jj=0; jj<cells_per_dim[1]; jj++) {
125c4762a1bSJed Brown       for (ii=0; ii<cells_per_dim[0]; ii++) {
126c4762a1bSJed Brown         PetscInt idx,idx0,idx1,idx2,idx3;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown         idx = (ii) + (jj) * nx[0];
129c4762a1bSJed Brown         idx0 = idx;
130c4762a1bSJed Brown         idx1 = idx0 + 1;
131c4762a1bSJed Brown         idx2 = idx1 + nx[0];
132c4762a1bSJed Brown         idx3 = idx0 + nx[0];
133c4762a1bSJed Brown 
134c4762a1bSJed Brown         tricells[3*cnt+0] = idx0;
135c4762a1bSJed Brown         tricells[3*cnt+1] = idx1;
136c4762a1bSJed Brown         tricells[3*cnt+2] = idx2;
137c4762a1bSJed Brown         cnt++;
138c4762a1bSJed Brown 
139c4762a1bSJed Brown         tricells[3*cnt+0] = idx0;
140c4762a1bSJed Brown         tricells[3*cnt+1] = idx2;
141c4762a1bSJed Brown         tricells[3*cnt+2] = idx3;
142c4762a1bSJed Brown         cnt++;
143c4762a1bSJed Brown       }
144c4762a1bSJed Brown     }
1459566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD,dim,n_tricells,n_trivert,3,PETSC_TRUE,tricells,dim,trivert,&celldm));
1469566063dSJacob Faibussowitsch     PetscCall(PetscFree(trivert));
1479566063dSJacob Faibussowitsch     PetscCall(PetscFree(tricells));
148c4762a1bSJed Brown   }
149*08401ef6SPierre Jolivet   PetscCheck(dim != 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only 2D PLEX example supported");
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Distribute mesh over processes */
1529566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh));
153c4762a1bSJed Brown   if (distributedMesh) {
1549566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&celldm));
155c4762a1bSJed Brown     celldm = distributedMesh;
156c4762a1bSJed Brown   }
1579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)celldm,"Cells"));
1589566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(celldm));
159c4762a1bSJed Brown   {
160c4762a1bSJed Brown     PetscInt     numComp[] = {1};
161c4762a1bSJed Brown     PetscInt     numDof[] = {1,0,0}; /* vert, edge, cell */
162c4762a1bSJed Brown     PetscInt     numBC = 0;
163c4762a1bSJed Brown     PetscSection section;
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section));
1669566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(celldm,section));
1679566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
168c4762a1bSJed Brown   }
1699566063dSJacob Faibussowitsch   PetscCall(DMSetUp(celldm));
170c4762a1bSJed Brown   {
171c4762a1bSJed Brown     PetscViewer viewer;
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
1749566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
1759566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
1769566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(viewer,"ex20plex.vtk"));
1779566063dSJacob Faibussowitsch     PetscCall(DMView(celldm,viewer));
1789566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Create the DMSwarm */
1829566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm));
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)swarm,"Swarm"));
1849566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm,DMSWARM));
1859566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm,dim));
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC));
1889566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm,celldm));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
1919566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE));
1929566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE));
1939566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
1969566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm,4,0));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,2));
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",1,fieldnames));
2019566063dSJacob Faibussowitsch   PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
2029566063dSJacob Faibussowitsch   PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
2039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
2049566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown PetscErrorCode pic_insert_DMPLEX(PetscBool is_simplex,PetscInt dim)
209c4762a1bSJed Brown {
210c4762a1bSJed Brown   DM             celldm,swarm,distributedMesh = NULL;
211c4762a1bSJed Brown   const char     *fieldnames[] = {"viscosity","DMSwarm_rank"};
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   PetscFunctionBegin;
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   /* Create the background cell DM */
216c4762a1bSJed Brown   {
217c4762a1bSJed Brown     PetscInt faces[3] = {4, 2, 4};
2189566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Distribute mesh over processes */
2229566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh));
223c4762a1bSJed Brown   if (distributedMesh) {
2249566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&celldm));
225c4762a1bSJed Brown     celldm = distributedMesh;
226c4762a1bSJed Brown   }
2279566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)celldm,"Cells"));
2289566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(celldm));
229c4762a1bSJed Brown   {
230c4762a1bSJed Brown     PetscInt     numComp[] = {1};
231c4762a1bSJed Brown     PetscInt     numDof[] = {1,0,0}; /* vert, edge, cell */
232c4762a1bSJed Brown     PetscInt     numBC = 0;
233c4762a1bSJed Brown     PetscSection section;
234c4762a1bSJed Brown 
2359566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section));
2369566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(celldm,section));
2379566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
238c4762a1bSJed Brown   }
2399566063dSJacob Faibussowitsch   PetscCall(DMSetUp(celldm));
240c4762a1bSJed Brown   {
241c4762a1bSJed Brown     PetscViewer viewer;
242c4762a1bSJed Brown 
2439566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
2459566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(viewer,"ex20plex.vtk"));
2479566063dSJacob Faibussowitsch     PetscCall(DMView(celldm,viewer));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
249c4762a1bSJed Brown   }
250c4762a1bSJed Brown 
2519566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm));
2529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)swarm,"Swarm"));
2539566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm,DMSWARM));
2549566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm,dim));
255c4762a1bSJed Brown 
2569566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC));
2579566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm,celldm));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
2609566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE));
2619566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE));
2629566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
2659566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm,4,0));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
2689566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_GAUSS,3));
2699566063dSJacob Faibussowitsch   PetscCall(DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",2,fieldnames));
2709566063dSJacob Faibussowitsch   PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
2719566063dSJacob Faibussowitsch   PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
2729566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
2739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
274c4762a1bSJed Brown   PetscFunctionReturn(0);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
277c4762a1bSJed Brown int main(int argc,char **args)
278c4762a1bSJed Brown {
279c4762a1bSJed Brown   PetscInt       mode = 0;
280c4762a1bSJed Brown   PetscInt       dim = 2;
281c4762a1bSJed Brown 
2829566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
2839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL));
2849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
285c4762a1bSJed Brown   switch (mode) {
286c4762a1bSJed Brown     case 0:
2879566063dSJacob Faibussowitsch       PetscCall(pic_insert_DMDA(dim));
288c4762a1bSJed Brown       break;
289c4762a1bSJed Brown     case 1:
290c4762a1bSJed Brown       /* tri / tet */
2919566063dSJacob Faibussowitsch       PetscCall(pic_insert_DMPLEX(PETSC_TRUE,dim));
292c4762a1bSJed Brown       break;
293c4762a1bSJed Brown     case 2:
294c4762a1bSJed Brown       /* quad / hex */
2959566063dSJacob Faibussowitsch       PetscCall(pic_insert_DMPLEX(PETSC_FALSE,dim));
296c4762a1bSJed Brown       break;
297c4762a1bSJed Brown     default:
2989566063dSJacob Faibussowitsch       PetscCall(pic_insert_DMDA(dim));
299c4762a1bSJed Brown       break;
300c4762a1bSJed Brown   }
3019566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
302b122ec5aSJacob Faibussowitsch   return 0;
303c4762a1bSJed Brown }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown /*TEST
306c4762a1bSJed Brown 
307c4762a1bSJed Brown    test:
308c4762a1bSJed Brown       args:
309c4762a1bSJed Brown       requires: !complex double
3107ed4f029SJed Brown       filter: grep -v atomic
311c4762a1bSJed Brown       filter_output: grep -v atomic
312c4762a1bSJed Brown 
313c4762a1bSJed Brown    test:
314c4762a1bSJed Brown       suffix: 2
315c4762a1bSJed Brown       requires: triangle double !complex
316c4762a1bSJed Brown       args: -mode 1
3177ed4f029SJed Brown       filter: grep -v atomic
318c4762a1bSJed Brown       filter_output: grep -v atomic
319c4762a1bSJed Brown 
320c4762a1bSJed Brown TEST*/
321