xref: /petsc/src/dm/tutorials/ex20.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "DMSwarm-PIC demonstator of inserting points into a cell DM \n\
3*c4762a1bSJed Brown Options: \n\
4*c4762a1bSJed Brown -mode {0,1} : 0 ==> DMDA, 1 ==> DMPLEX cell DM \n\
5*c4762a1bSJed Brown -dim {2,3}  : spatial dimension\n";
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown #include <petsc.h>
8*c4762a1bSJed Brown #include <petscdm.h>
9*c4762a1bSJed Brown #include <petscdmda.h>
10*c4762a1bSJed Brown #include <petscdmplex.h>
11*c4762a1bSJed Brown #include <petscdmswarm.h>
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown PetscErrorCode pic_insert_DMDA(PetscInt dim)
14*c4762a1bSJed Brown {
15*c4762a1bSJed Brown   PetscErrorCode ierr;
16*c4762a1bSJed Brown   DM             celldm = NULL,swarm;
17*c4762a1bSJed Brown   PetscInt       dof,stencil_width;
18*c4762a1bSJed Brown   PetscReal      min[3],max[3];
19*c4762a1bSJed Brown   PetscInt       ndir[3];
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   PetscFunctionBegin;
22*c4762a1bSJed Brown   /* Create the background cell DM */
23*c4762a1bSJed Brown   dof = 1;
24*c4762a1bSJed Brown   stencil_width = 1;
25*c4762a1bSJed Brown   if (dim == 2) {
26*c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
27*c4762a1bSJed Brown   }
28*c4762a1bSJed Brown   if (dim == 3) {
29*c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
30*c4762a1bSJed Brown   }
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown   ierr = DMDASetElementType(celldm,DMDA_ELEMENT_Q1);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = DMSetFromOptions(celldm);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = DMSetUp(celldm);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(celldm,0.0,2.0,0.0,1.0,0.0,1.5);CHKERRQ(ierr);
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   /* Create the DMSwarm */
39*c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&swarm);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = DMSetType(swarm,DMSWARM);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = DMSetDimension(swarm,dim);CHKERRQ(ierr);
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
44*c4762a1bSJed Brown   ierr = DMSwarmSetType(swarm,DMSWARM_PIC);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(swarm,celldm);CHKERRQ(ierr);
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
48*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(swarm);CHKERRQ(ierr);
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
53*c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(swarm,4,0);CHKERRQ(ierr);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
56*c4762a1bSJed Brown   ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,3);CHKERRQ(ierr);
57*c4762a1bSJed Brown   min[0] = 0.5; max[0] = 0.7;
58*c4762a1bSJed Brown   min[1] = 0.5; max[1] = 0.8;
59*c4762a1bSJed Brown   min[2] = 0.5; max[2] = 0.9;
60*c4762a1bSJed Brown   ndir[0] = ndir[1] = ndir[2] = 30;
61*c4762a1bSJed Brown   ierr = DMSwarmSetPointsUniformCoordinates(swarm,min,max,ndir,ADD_VALUES);CHKERRQ(ierr);
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   /* This should be dispatched from a regular DMView() */
64*c4762a1bSJed Brown   ierr = DMSwarmViewXDMF(swarm,"ex20.xmf");CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = DMView(celldm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = DMView(swarm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   {
69*c4762a1bSJed Brown     PetscInt    npoints,*list;
70*c4762a1bSJed Brown     PetscMPIInt rank;
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown     ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
73*c4762a1bSJed Brown     ierr = DMSwarmSortGetAccess(swarm);CHKERRQ(ierr);
74*c4762a1bSJed Brown     ierr = DMSwarmSortGetNumberOfPointsPerCell(swarm,0,&npoints);CHKERRQ(ierr);
75*c4762a1bSJed Brown     ierr = DMSwarmSortGetPointsPerCell(swarm,rank,&npoints,&list);CHKERRQ(ierr);
76*c4762a1bSJed Brown     ierr = PetscFree(list);CHKERRQ(ierr);
77*c4762a1bSJed Brown     ierr = DMSwarmSortRestoreAccess(swarm);CHKERRQ(ierr);
78*c4762a1bSJed Brown   }
79*c4762a1bSJed Brown   ierr = DMSwarmMigrate(swarm,PETSC_FALSE);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = DMDestroy(&celldm);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = DMDestroy(&swarm);CHKERRQ(ierr);
82*c4762a1bSJed Brown   PetscFunctionReturn(0);
83*c4762a1bSJed Brown }
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown PetscErrorCode pic_insert_DMPLEX_with_cell_list(PetscInt dim)
86*c4762a1bSJed Brown {
87*c4762a1bSJed Brown   PetscErrorCode ierr;
88*c4762a1bSJed Brown   DM             celldm = NULL,swarm,distributedMesh = NULL;
89*c4762a1bSJed Brown   const  char    *fieldnames[] = {"viscosity"};
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown   PetscFunctionBegin;
92*c4762a1bSJed Brown   /* Create the background cell DM */
93*c4762a1bSJed Brown   if (dim == 2) {
94*c4762a1bSJed Brown     PetscInt cells_per_dim[2],nx[2];
95*c4762a1bSJed Brown     PetscInt n_tricells;
96*c4762a1bSJed Brown     PetscInt n_trivert;
97*c4762a1bSJed Brown     int      *tricells;
98*c4762a1bSJed Brown     double   *trivert,dx,dy;
99*c4762a1bSJed Brown     PetscInt ii,jj,cnt;
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown     cells_per_dim[0] = 4;
102*c4762a1bSJed Brown     cells_per_dim[1] = 4;
103*c4762a1bSJed Brown     n_tricells = cells_per_dim[0] * cells_per_dim[1] * 2;
104*c4762a1bSJed Brown     nx[0] = cells_per_dim[0] + 1;
105*c4762a1bSJed Brown     nx[1] = cells_per_dim[1] + 1;
106*c4762a1bSJed Brown     n_trivert = nx[0] * nx[1];
107*c4762a1bSJed Brown 
108*c4762a1bSJed Brown     ierr = PetscMalloc1(n_tricells*3,&tricells);CHKERRQ(ierr);
109*c4762a1bSJed Brown     ierr = PetscMalloc1(nx[0]*nx[1]*2,&trivert);CHKERRQ(ierr);
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown     /* verts */
112*c4762a1bSJed Brown     cnt = 0;
113*c4762a1bSJed Brown     dx = 2.0/((PetscReal)cells_per_dim[0]);
114*c4762a1bSJed Brown     dy = 1.0/((PetscReal)cells_per_dim[1]);
115*c4762a1bSJed Brown     for (jj=0; jj<nx[1]; jj++) {
116*c4762a1bSJed Brown       for (ii=0; ii<nx[0]; ii++) {
117*c4762a1bSJed Brown         trivert[2*cnt+0] = 0.0 + ii * dx;
118*c4762a1bSJed Brown         trivert[2*cnt+1] = 0.0 + jj * dy;
119*c4762a1bSJed Brown         cnt++;
120*c4762a1bSJed Brown       }
121*c4762a1bSJed Brown     }
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown     /* connectivity */
124*c4762a1bSJed Brown     cnt = 0;
125*c4762a1bSJed Brown     for (jj=0; jj<cells_per_dim[1]; jj++) {
126*c4762a1bSJed Brown       for (ii=0; ii<cells_per_dim[0]; ii++) {
127*c4762a1bSJed Brown         PetscInt idx,idx0,idx1,idx2,idx3;
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown         idx = (ii) + (jj) * nx[0];
130*c4762a1bSJed Brown         idx0 = idx;
131*c4762a1bSJed Brown         idx1 = idx0 + 1;
132*c4762a1bSJed Brown         idx2 = idx1 + nx[0];
133*c4762a1bSJed Brown         idx3 = idx0 + nx[0];
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown         tricells[3*cnt+0] = idx0;
136*c4762a1bSJed Brown         tricells[3*cnt+1] = idx1;
137*c4762a1bSJed Brown         tricells[3*cnt+2] = idx2;
138*c4762a1bSJed Brown         cnt++;
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown         tricells[3*cnt+0] = idx0;
141*c4762a1bSJed Brown         tricells[3*cnt+1] = idx2;
142*c4762a1bSJed Brown         tricells[3*cnt+2] = idx3;
143*c4762a1bSJed Brown         cnt++;
144*c4762a1bSJed Brown       }
145*c4762a1bSJed Brown     }
146*c4762a1bSJed Brown     ierr = DMPlexCreateFromCellList(PETSC_COMM_WORLD,dim,n_tricells,n_trivert,3,PETSC_TRUE,tricells,dim,trivert,&celldm);CHKERRQ(ierr);
147*c4762a1bSJed Brown     ierr = PetscFree(trivert);CHKERRQ(ierr);
148*c4762a1bSJed Brown     ierr = PetscFree(tricells);CHKERRQ(ierr);
149*c4762a1bSJed Brown   }
150*c4762a1bSJed Brown   if (dim == 3) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only 2D PLEX example supported");
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   /* Distribute mesh over processes */
153*c4762a1bSJed Brown   ierr = DMPlexDistribute(celldm,0,NULL,&distributedMesh);CHKERRQ(ierr);
154*c4762a1bSJed Brown   if (distributedMesh) {
155*c4762a1bSJed Brown     ierr = DMDestroy(&celldm);CHKERRQ(ierr);
156*c4762a1bSJed Brown     celldm = distributedMesh;
157*c4762a1bSJed Brown   }
158*c4762a1bSJed Brown   ierr = DMSetFromOptions(celldm);CHKERRQ(ierr);
159*c4762a1bSJed Brown   {
160*c4762a1bSJed Brown     PetscInt     numComp[] = {1};
161*c4762a1bSJed Brown     PetscInt     numDof[] = {1,0,0}; /* vert, edge, cell */
162*c4762a1bSJed Brown     PetscInt     numBC = 0;
163*c4762a1bSJed Brown     PetscSection section;
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown     ierr = DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section);CHKERRQ(ierr);
166*c4762a1bSJed Brown     ierr = DMSetLocalSection(celldm,section);CHKERRQ(ierr);
167*c4762a1bSJed Brown     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
168*c4762a1bSJed Brown   }
169*c4762a1bSJed Brown   ierr = DMSetUp(celldm);CHKERRQ(ierr);
170*c4762a1bSJed Brown   {
171*c4762a1bSJed Brown     PetscViewer viewer;
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
174*c4762a1bSJed Brown     ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
175*c4762a1bSJed Brown     ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
176*c4762a1bSJed Brown     ierr = PetscViewerFileSetName(viewer,"ex20plex.vtk");CHKERRQ(ierr);
177*c4762a1bSJed Brown     ierr = DMView(celldm,viewer);CHKERRQ(ierr);
178*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
179*c4762a1bSJed Brown   }
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   /* Create the DMSwarm */
182*c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&swarm);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = DMSetType(swarm,DMSWARM);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = DMSetDimension(swarm,dim);CHKERRQ(ierr);
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown   ierr = DMSwarmSetType(swarm,DMSWARM_PIC);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(swarm,celldm);CHKERRQ(ierr);
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
190*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(swarm);CHKERRQ(ierr);
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
195*c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(swarm,4,0);CHKERRQ(ierr);
196*c4762a1bSJed Brown 
197*c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
198*c4762a1bSJed Brown   ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,2);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr = DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",1,fieldnames);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = DMView(celldm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = DMView(swarm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = DMDestroy(&celldm);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = DMDestroy(&swarm);CHKERRQ(ierr);
204*c4762a1bSJed Brown   PetscFunctionReturn(0);
205*c4762a1bSJed Brown }
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown PetscErrorCode pic_insert_DMPLEX(PetscBool is_simplex,PetscInt dim)
208*c4762a1bSJed Brown {
209*c4762a1bSJed Brown   PetscErrorCode ierr;
210*c4762a1bSJed Brown   DM             celldm,swarm,distributedMesh = NULL;
211*c4762a1bSJed Brown   const char     *fieldnames[] = {"viscosity","DMSwarm_rank"};
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown   PetscFunctionBegin;
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   /* Create the background cell DM */
216*c4762a1bSJed Brown   {
217*c4762a1bSJed Brown     PetscInt faces[3] = {4, 2, 4};
218*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm);CHKERRQ(ierr);
219*c4762a1bSJed Brown   }
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   /* Distribute mesh over processes */
222*c4762a1bSJed Brown   ierr = DMPlexDistribute(celldm,0,NULL,&distributedMesh);CHKERRQ(ierr);
223*c4762a1bSJed Brown   if (distributedMesh) {
224*c4762a1bSJed Brown     ierr = DMDestroy(&celldm);CHKERRQ(ierr);
225*c4762a1bSJed Brown     celldm = distributedMesh;
226*c4762a1bSJed Brown   }
227*c4762a1bSJed Brown   ierr = DMSetFromOptions(celldm);CHKERRQ(ierr);
228*c4762a1bSJed Brown   {
229*c4762a1bSJed Brown     PetscInt     numComp[] = {1};
230*c4762a1bSJed Brown     PetscInt     numDof[] = {1,0,0}; /* vert, edge, cell */
231*c4762a1bSJed Brown     PetscInt     numBC = 0;
232*c4762a1bSJed Brown     PetscSection section;
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown     ierr = DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section);CHKERRQ(ierr);
235*c4762a1bSJed Brown     ierr = DMSetLocalSection(celldm,section);CHKERRQ(ierr);
236*c4762a1bSJed Brown     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
237*c4762a1bSJed Brown   }
238*c4762a1bSJed Brown   ierr = DMSetUp(celldm);CHKERRQ(ierr);
239*c4762a1bSJed Brown   {
240*c4762a1bSJed Brown     PetscViewer viewer;
241*c4762a1bSJed Brown 
242*c4762a1bSJed Brown     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
243*c4762a1bSJed Brown     ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
244*c4762a1bSJed Brown     ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
245*c4762a1bSJed Brown     ierr = PetscViewerFileSetName(viewer,"ex20plex.vtk");CHKERRQ(ierr);
246*c4762a1bSJed Brown     ierr = DMView(celldm,viewer);CHKERRQ(ierr);
247*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
248*c4762a1bSJed Brown   }
249*c4762a1bSJed Brown 
250*c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&swarm);CHKERRQ(ierr);
251*c4762a1bSJed Brown   ierr = DMSetType(swarm,DMSWARM);CHKERRQ(ierr);
252*c4762a1bSJed Brown   ierr = DMSetDimension(swarm,dim);CHKERRQ(ierr);
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   ierr = DMSwarmSetType(swarm,DMSWARM_PIC);CHKERRQ(ierr);
255*c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(swarm,celldm);CHKERRQ(ierr);
256*c4762a1bSJed Brown 
257*c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
258*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE);CHKERRQ(ierr);
259*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE);CHKERRQ(ierr);
260*c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(swarm);CHKERRQ(ierr);
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
263*c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(swarm,4,0);CHKERRQ(ierr);
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
266*c4762a1bSJed Brown   ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_GAUSS,3);CHKERRQ(ierr);
267*c4762a1bSJed Brown   ierr = DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",2,fieldnames);CHKERRQ(ierr);
268*c4762a1bSJed Brown   ierr = DMView(celldm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = DMView(swarm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
270*c4762a1bSJed Brown   ierr = DMDestroy(&celldm);CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = DMDestroy(&swarm);CHKERRQ(ierr);
272*c4762a1bSJed Brown   PetscFunctionReturn(0);
273*c4762a1bSJed Brown }
274*c4762a1bSJed Brown 
275*c4762a1bSJed Brown int main(int argc,char **args)
276*c4762a1bSJed Brown {
277*c4762a1bSJed Brown   PetscErrorCode ierr;
278*c4762a1bSJed Brown   PetscInt       mode = 0;
279*c4762a1bSJed Brown   PetscInt       dim = 2;
280*c4762a1bSJed Brown 
281*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
282*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL);CHKERRQ(ierr);
283*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr);
284*c4762a1bSJed Brown   switch (mode) {
285*c4762a1bSJed Brown     case 0:
286*c4762a1bSJed Brown       ierr = pic_insert_DMDA(dim);CHKERRQ(ierr);
287*c4762a1bSJed Brown       break;
288*c4762a1bSJed Brown     case 1:
289*c4762a1bSJed Brown       /* tri / tet */
290*c4762a1bSJed Brown       ierr = pic_insert_DMPLEX(PETSC_TRUE,dim);CHKERRQ(ierr);
291*c4762a1bSJed Brown       break;
292*c4762a1bSJed Brown     case 2:
293*c4762a1bSJed Brown       /* quad / hex */
294*c4762a1bSJed Brown       ierr = pic_insert_DMPLEX(PETSC_FALSE,dim);CHKERRQ(ierr);
295*c4762a1bSJed Brown       break;
296*c4762a1bSJed Brown     default:
297*c4762a1bSJed Brown       ierr = pic_insert_DMDA(dim);CHKERRQ(ierr);
298*c4762a1bSJed Brown       break;
299*c4762a1bSJed Brown   }
300*c4762a1bSJed Brown   ierr = PetscFinalize();
301*c4762a1bSJed Brown   return ierr;
302*c4762a1bSJed Brown }
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown /*TEST
305*c4762a1bSJed Brown 
306*c4762a1bSJed Brown    test:
307*c4762a1bSJed Brown       args:
308*c4762a1bSJed Brown       requires: !complex double
309*c4762a1bSJed Brown       filter: grep -v DM_ | grep -v atomic
310*c4762a1bSJed Brown       filter_output: grep -v atomic
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown    test:
313*c4762a1bSJed Brown       suffix: 2
314*c4762a1bSJed Brown       requires: triangle double !complex
315*c4762a1bSJed Brown       args: -mode 1
316*c4762a1bSJed Brown       filter: grep -v DM_ | grep -v atomic
317*c4762a1bSJed Brown       filter_output: grep -v atomic
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown TEST*/
320