xref: /petsc/src/dm/tutorials/ex20.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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) {
25*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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) {
28*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetElementType(celldm,DMDA_ELEMENT_Q1));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(celldm));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(celldm));
34c4762a1bSJed Brown 
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(celldm,0.0,2.0,0.0,1.0,0.0,1.5));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Create the DMSwarm */
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD,&swarm));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)swarm,"Swarm"));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(swarm,DMSWARM));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(swarm,dim));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(swarm,DMSWARM_PIC));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(swarm,celldm));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(swarm));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(swarm,4,0));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetPointsUniformCoordinates(swarm,min,max,ndir,ADD_VALUES));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* This should be dispatched from a regular DMView() */
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmViewXDMF(swarm,"ex20.xmf"));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   {
69c4762a1bSJed Brown     PetscInt    npoints,*list;
70c4762a1bSJed Brown     PetscMPIInt rank;
71c4762a1bSJed Brown 
72*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortGetAccess(swarm));
74*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortGetNumberOfPointsPerCell(swarm,0,&npoints));
75*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortGetPointsPerCell(swarm,rank,&npoints,&list));
76*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(list));
77*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortRestoreAccess(swarm));
78c4762a1bSJed Brown   }
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmMigrate(swarm,PETSC_FALSE));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&celldm));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
107*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(n_tricells*3,&tricells));
108*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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     }
145*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD,dim,n_tricells,n_trivert,3,PETSC_TRUE,tricells,dim,trivert,&celldm));
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(trivert));
147*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(tricells));
148c4762a1bSJed Brown   }
1492c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim == 3,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only 2D PLEX example supported");
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Distribute mesh over processes */
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(celldm,0,NULL,&distributedMesh));
153c4762a1bSJed Brown   if (distributedMesh) {
154*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&celldm));
155c4762a1bSJed Brown     celldm = distributedMesh;
156c4762a1bSJed Brown   }
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)celldm,"Cells"));
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
165*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section));
166*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetLocalSection(celldm,section));
167*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionDestroy(&section));
168c4762a1bSJed Brown   }
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(celldm));
170c4762a1bSJed Brown   {
171c4762a1bSJed Brown     PetscViewer viewer;
172c4762a1bSJed Brown 
173*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
174*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWERVTK));
175*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
176*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFileSetName(viewer,"ex20plex.vtk"));
177*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMView(celldm,viewer));
178*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* Create the DMSwarm */
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD,&swarm));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)swarm,"Swarm"));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(swarm,DMSWARM));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(swarm,dim));
186c4762a1bSJed Brown 
187*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(swarm,DMSWARM_PIC));
188*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(swarm,celldm));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE));
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE));
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(swarm));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(swarm,4,0));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,2));
200*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",1,fieldnames));
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&celldm));
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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};
218*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, is_simplex, faces, NULL, NULL, NULL, PETSC_TRUE, &celldm));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* Distribute mesh over processes */
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexDistribute(celldm,0,NULL,&distributedMesh));
223c4762a1bSJed Brown   if (distributedMesh) {
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&celldm));
225c4762a1bSJed Brown     celldm = distributedMesh;
226c4762a1bSJed Brown   }
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)celldm,"Cells"));
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
235*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section));
236*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetLocalSection(celldm,section));
237*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionDestroy(&section));
238c4762a1bSJed Brown   }
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(celldm));
240c4762a1bSJed Brown   {
241c4762a1bSJed Brown     PetscViewer viewer;
242c4762a1bSJed Brown 
243*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
244*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWERVTK));
245*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
246*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerFileSetName(viewer,"ex20plex.vtk"));
247*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMView(celldm,viewer));
248*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
249c4762a1bSJed Brown   }
250c4762a1bSJed Brown 
251*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD,&swarm));
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)swarm,"Swarm"));
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(swarm,DMSWARM));
254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(swarm,dim));
255c4762a1bSJed Brown 
256*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(swarm,DMSWARM_PIC));
257*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(swarm,celldm));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
260*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"viscosity",1,PETSC_DOUBLE));
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(swarm,"density",1,PETSC_DOUBLE));
262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(swarm));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(swarm,4,0));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_GAUSS,3));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmViewFieldsXDMF(swarm,"ex20.xmf",2,fieldnames));
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
272*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&celldm));
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&swarm));
274c4762a1bSJed Brown   PetscFunctionReturn(0);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
277c4762a1bSJed Brown int main(int argc,char **args)
278c4762a1bSJed Brown {
279c4762a1bSJed Brown   PetscErrorCode ierr;
280c4762a1bSJed Brown   PetscInt       mode = 0;
281c4762a1bSJed Brown   PetscInt       dim = 2;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL));
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
286c4762a1bSJed Brown   switch (mode) {
287c4762a1bSJed Brown     case 0:
288*5f80ce2aSJacob Faibussowitsch       CHKERRQ(pic_insert_DMDA(dim));
289c4762a1bSJed Brown       break;
290c4762a1bSJed Brown     case 1:
291c4762a1bSJed Brown       /* tri / tet */
292*5f80ce2aSJacob Faibussowitsch       CHKERRQ(pic_insert_DMPLEX(PETSC_TRUE,dim));
293c4762a1bSJed Brown       break;
294c4762a1bSJed Brown     case 2:
295c4762a1bSJed Brown       /* quad / hex */
296*5f80ce2aSJacob Faibussowitsch       CHKERRQ(pic_insert_DMPLEX(PETSC_FALSE,dim));
297c4762a1bSJed Brown       break;
298c4762a1bSJed Brown     default:
299*5f80ce2aSJacob Faibussowitsch       CHKERRQ(pic_insert_DMDA(dim));
300c4762a1bSJed Brown       break;
301c4762a1bSJed Brown   }
302c4762a1bSJed Brown   ierr = PetscFinalize();
303c4762a1bSJed Brown   return ierr;
304c4762a1bSJed Brown }
305c4762a1bSJed Brown 
306c4762a1bSJed Brown /*TEST
307c4762a1bSJed Brown 
308c4762a1bSJed Brown    test:
309c4762a1bSJed Brown       args:
310c4762a1bSJed Brown       requires: !complex double
3117ed4f029SJed Brown       filter: grep -v atomic
312c4762a1bSJed Brown       filter_output: grep -v atomic
313c4762a1bSJed Brown 
314c4762a1bSJed Brown    test:
315c4762a1bSJed Brown       suffix: 2
316c4762a1bSJed Brown       requires: triangle double !complex
317c4762a1bSJed Brown       args: -mode 1
3187ed4f029SJed Brown       filter: grep -v atomic
319c4762a1bSJed Brown       filter_output: grep -v atomic
320c4762a1bSJed Brown 
321c4762a1bSJed Brown TEST*/
322