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