xref: /petsc/src/dm/tutorials/ex21.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "DMSwarm-PIC demonstator of advecting points within cell DM defined by a DA or PLEX object \n\
3c4762a1bSJed Brown Options: \n\
4c4762a1bSJed Brown -ppcell   : Number of times to sub-divide the reference cell when layout the initial particle coordinates \n\
5c4762a1bSJed Brown -meshtype : 0 ==> DA , 1 ==> PLEX \n\
6c4762a1bSJed Brown -nt       : Number of timestep to perform \n\
7c4762a1bSJed Brown -view     : Write out initial condition and time dependent data \n";
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include <petsc.h>
10c4762a1bSJed Brown #include <petscdm.h>
11c4762a1bSJed Brown #include <petscdmda.h>
12c4762a1bSJed Brown #include <petscdmswarm.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown PetscErrorCode pic_advect(PetscInt ppcell,PetscInt meshtype)
15c4762a1bSJed Brown {
16c4762a1bSJed Brown   PetscErrorCode ierr;
17c4762a1bSJed Brown   const PetscInt dim = 2;
18c4762a1bSJed Brown   DM celldm,swarm;
19c4762a1bSJed Brown   PetscInt tk,nt = 200;
20c4762a1bSJed Brown   PetscBool view = PETSC_FALSE;
21c4762a1bSJed Brown   Vec *pfields;
22c4762a1bSJed Brown   PetscReal minradius;
23c4762a1bSJed Brown   PetscReal dt;
24c4762a1bSJed Brown   PetscReal vel[] = { 1.0, 0.16 };
25c4762a1bSJed Brown   const char *fieldnames[] = { "phi" };
26c4762a1bSJed Brown   PetscViewer viewer;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   PetscFunctionBegin;
29c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nt",&nt,NULL);CHKERRQ(ierr);
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Create the background cell DM */
33c4762a1bSJed Brown   if (meshtype == 0) { /* DA */
34c4762a1bSJed Brown     PetscInt nxy;
35c4762a1bSJed Brown     PetscInt dof = 1;
36c4762a1bSJed Brown     PetscInt stencil_width = 1;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMDA\n");CHKERRQ(ierr);
39c4762a1bSJed Brown     nxy = 33;
40c4762a1bSJed Brown     ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,nxy,nxy,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&celldm);CHKERRQ(ierr);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown     ierr = DMDASetElementType(celldm,DMDA_ELEMENT_Q1);CHKERRQ(ierr);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown     ierr = DMSetFromOptions(celldm);CHKERRQ(ierr);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown     ierr = DMSetUp(celldm);CHKERRQ(ierr);
47c4762a1bSJed Brown 
48c4762a1bSJed Brown     ierr = DMDASetUniformCoordinates(celldm,0.0,1.0,0.0,1.0,0.0,1.5);CHKERRQ(ierr);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown     minradius = 1.0/((PetscReal)(nxy-1));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"DA(minradius) %1.4e\n",(double)minradius);CHKERRQ(ierr);
53c4762a1bSJed Brown   }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   if (meshtype == 1){ /* PLEX */
56c4762a1bSJed Brown     DM distributedMesh = NULL;
57c4762a1bSJed Brown     PetscInt numComp[] = {1};
58c4762a1bSJed Brown     PetscInt numDof[] = {1,0,0}; /* vert, edge, cell */
59c4762a1bSJed Brown     PetscInt faces[]  = {1,1,1};
60c4762a1bSJed Brown     PetscInt numBC = 0;
61c4762a1bSJed Brown     PetscSection section;
62c4762a1bSJed Brown     Vec cellgeom = NULL;
63c4762a1bSJed Brown     Vec facegeom = NULL;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMPLEX\n");CHKERRQ(ierr);
66c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, &celldm);CHKERRQ(ierr);
67c4762a1bSJed Brown 
68c4762a1bSJed Brown     /* Distribute mesh over processes */
69c4762a1bSJed Brown     ierr = DMPlexDistribute(celldm,0,NULL,&distributedMesh);CHKERRQ(ierr);
70c4762a1bSJed Brown     if (distributedMesh) {
71c4762a1bSJed Brown       ierr = DMDestroy(&celldm);CHKERRQ(ierr);
72c4762a1bSJed Brown       celldm = distributedMesh;
73c4762a1bSJed Brown     }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     ierr = DMSetFromOptions(celldm);CHKERRQ(ierr);
76c4762a1bSJed Brown 
77c4762a1bSJed Brown     ierr = DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section);CHKERRQ(ierr);
78c4762a1bSJed Brown     ierr = DMSetLocalSection(celldm,section);CHKERRQ(ierr);
79c4762a1bSJed Brown 
80c4762a1bSJed Brown     ierr = DMSetUp(celldm);CHKERRQ(ierr);
81c4762a1bSJed Brown 
82c4762a1bSJed Brown     /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */
83c4762a1bSJed Brown     ierr = DMPlexComputeGeometryFVM(celldm,&cellgeom,&facegeom);CHKERRQ(ierr);
84c4762a1bSJed Brown     ierr = DMPlexGetMinRadius(celldm,&minradius);CHKERRQ(ierr);
85c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"PLEX(minradius) %1.4e\n",(double)minradius);CHKERRQ(ierr);
86c4762a1bSJed Brown     ierr = VecDestroy(&cellgeom);CHKERRQ(ierr);
87c4762a1bSJed Brown     ierr = VecDestroy(&facegeom);CHKERRQ(ierr);
88c4762a1bSJed Brown     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Create the DMSwarm */
92c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&swarm);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = DMSetType(swarm,DMSWARM);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = DMSetDimension(swarm,dim);CHKERRQ(ierr);
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
97c4762a1bSJed Brown   ierr = DMSwarmSetType(swarm,DMSWARM_PIC);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(swarm,celldm);CHKERRQ(ierr);
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
101c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"phi",1,PETSC_REAL);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(swarm,"region",1,PETSC_REAL);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(swarm);CHKERRQ(ierr);
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
106c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(swarm,4,0);CHKERRQ(ierr);
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
109c4762a1bSJed Brown   /*ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell);CHKERRQ(ierr);*/
110c4762a1bSJed Brown   ierr = DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,ppcell);CHKERRQ(ierr);
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* Define initial conditions for th swarm fields "phi" and "region" */
113c4762a1bSJed Brown   {
114c4762a1bSJed Brown     PetscReal *s_coor,*s_phi,*s_region;
115c4762a1bSJed Brown     PetscInt npoints,p;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr);
118c4762a1bSJed Brown     ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor);CHKERRQ(ierr);
119c4762a1bSJed Brown     ierr = DMSwarmGetField(swarm,"phi",NULL,NULL,(void**)&s_phi);CHKERRQ(ierr);
120c4762a1bSJed Brown     ierr = DMSwarmGetField(swarm,"region",NULL,NULL,(void**)&s_region);CHKERRQ(ierr);
121c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
122c4762a1bSJed Brown       PetscReal pos[2];
123c4762a1bSJed Brown       pos[0] = s_coor[2*p+0];
124c4762a1bSJed Brown       pos[1] = s_coor[2*p+1];
125c4762a1bSJed Brown 
126c4762a1bSJed Brown       s_region[p] = 1.0;
127c4762a1bSJed Brown       s_phi[p] = 1.0 + PetscExpReal(-200.0*((pos[0]-0.5)*(pos[0]-0.5) + (pos[1]-0.5)*(pos[1]-0.5)));
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown     ierr = DMSwarmRestoreField(swarm,"region",NULL,NULL,(void**)&s_region);CHKERRQ(ierr);
130c4762a1bSJed Brown     ierr = DMSwarmRestoreField(swarm,"phi",NULL,NULL,(void**)&s_phi);CHKERRQ(ierr);
131c4762a1bSJed Brown     ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor);CHKERRQ(ierr);
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Project initial value of phi onto the mesh */
135c4762a1bSJed Brown   ierr = DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_FALSE);CHKERRQ(ierr);
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   if (view) {
138c4762a1bSJed Brown     /* View swarm all swarm fields using data type PETSC_REAL */
139c4762a1bSJed Brown     ierr = DMSwarmViewXDMF(swarm,"ic_dms.xmf");CHKERRQ(ierr);
140c4762a1bSJed Brown 
141c4762a1bSJed Brown     /* View projected swarm field "phi" */
142c4762a1bSJed Brown     ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
143c4762a1bSJed Brown     ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
144c4762a1bSJed Brown     ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
145c4762a1bSJed Brown     if (meshtype == 0) { /* DA */
146c4762a1bSJed Brown       ierr = PetscViewerFileSetName(viewer,"ic_dmda.vts");CHKERRQ(ierr);
147c4762a1bSJed Brown       ierr = VecView(pfields[0],viewer);CHKERRQ(ierr);
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown     if (meshtype == 1) { /* PLEX */
150c4762a1bSJed Brown       ierr = PetscViewerFileSetName(viewer,"ic_dmplex.vtk");CHKERRQ(ierr);
151c4762a1bSJed Brown       ierr = DMView(celldm,viewer);CHKERRQ(ierr);
152c4762a1bSJed Brown       ierr = VecView(pfields[0],viewer);CHKERRQ(ierr);
153c4762a1bSJed Brown     }
154c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   ierr = DMView(celldm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = DMView(swarm,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   dt = 0.5 * minradius / PetscSqrtReal(vel[0]*vel[0] + vel[1]*vel[1]);
161c4762a1bSJed Brown   for (tk=1; tk<=nt; tk++) {
162c4762a1bSJed Brown     PetscReal *s_coor;
163c4762a1bSJed Brown     PetscInt npoints,p;
164c4762a1bSJed Brown 
165b5675b0fSBarry Smith     ierr = PetscPrintf(PETSC_COMM_WORLD,"[step %D]\n",tk);CHKERRQ(ierr);
166c4762a1bSJed Brown     /* advect with analytic prescribed (constant) velocity field */
167c4762a1bSJed Brown     ierr = DMSwarmGetLocalSize(swarm,&npoints);CHKERRQ(ierr);
168c4762a1bSJed Brown     ierr = DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor);CHKERRQ(ierr);
169c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
170c4762a1bSJed Brown       s_coor[2*p+0] += dt * vel[0];
171c4762a1bSJed Brown       s_coor[2*p+1] += dt * vel[1];
172c4762a1bSJed Brown     }
173c4762a1bSJed Brown     ierr = DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor);CHKERRQ(ierr);
174c4762a1bSJed Brown 
175c4762a1bSJed Brown     ierr = DMSwarmMigrate(swarm,PETSC_TRUE);CHKERRQ(ierr);
176c4762a1bSJed Brown 
177c4762a1bSJed Brown     /* Ad-hoc cell filling algorithm */
178c4762a1bSJed Brown     /*
179c4762a1bSJed Brown        The injection frequency is chosen for default DA case.
180c4762a1bSJed Brown        They will likely not be appropriate for the general case using an unstructure PLEX mesh.
181c4762a1bSJed Brown     */
182c4762a1bSJed Brown     if (tk%10 == 0) {
183c4762a1bSJed Brown       PetscReal dx = 1.0/32.0;
184c4762a1bSJed Brown       PetscInt npoints_dir_x[] = { 32, 1 };
185c4762a1bSJed Brown       PetscReal min[2],max[2];
186c4762a1bSJed Brown 
187c4762a1bSJed Brown       min[0] = 0.5 * dx;  max[0] = 0.5 * dx + 31.0 * dx;
188c4762a1bSJed Brown       min[1] = 0.5 * dx;  max[1] = 0.5 * dx;
189c4762a1bSJed Brown       ierr = DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_x,ADD_VALUES);CHKERRQ(ierr);
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown     if (tk%2 == 0) {
192c4762a1bSJed Brown       PetscReal dx = 1.0/32.0;
193c4762a1bSJed Brown       PetscInt npoints_dir_y[] = { 2, 31 };
194c4762a1bSJed Brown       PetscReal min[2],max[2];
195c4762a1bSJed Brown 
196c4762a1bSJed Brown       min[0] = 0.05 * dx; max[0] = 0.5 * dx;
197c4762a1bSJed Brown       min[1] = 0.5 * dx;  max[1] = 0.5 * dx + 31.0 * dx;
198c4762a1bSJed Brown       ierr = DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_y,ADD_VALUES);CHKERRQ(ierr);
199c4762a1bSJed Brown     }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown     /* Project swarm field "phi" onto the cell DM */
202c4762a1bSJed Brown     ierr = DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_TRUE);CHKERRQ(ierr);
203c4762a1bSJed Brown 
204c4762a1bSJed Brown     if (view) {
205c4762a1bSJed Brown       PetscViewer viewer;
206c4762a1bSJed Brown       char fname[PETSC_MAX_PATH_LEN];
207c4762a1bSJed Brown 
208c4762a1bSJed Brown       /* View swarm fields */
209c4762a1bSJed Brown       ierr = PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dms.xmf",tk);CHKERRQ(ierr);
210c4762a1bSJed Brown       ierr = DMSwarmViewXDMF(swarm,fname);CHKERRQ(ierr);
211c4762a1bSJed Brown 
212c4762a1bSJed Brown       /* View projected field */
213c4762a1bSJed Brown       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
214c4762a1bSJed Brown       ierr = PetscViewerSetType(viewer,PETSCVIEWERVTK);CHKERRQ(ierr);
215c4762a1bSJed Brown       ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
216c4762a1bSJed Brown 
217c4762a1bSJed Brown       if (meshtype == 0) { /* DA */
218c4762a1bSJed Brown         ierr = PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmda.vts",tk);CHKERRQ(ierr);
219c4762a1bSJed Brown         ierr = PetscViewerFileSetName(viewer,fname);CHKERRQ(ierr);
220c4762a1bSJed Brown         ierr = VecView(pfields[0],viewer);CHKERRQ(ierr);
221c4762a1bSJed Brown       }
222c4762a1bSJed Brown       if (meshtype == 1) { /* PLEX */
223c4762a1bSJed Brown         ierr = PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmplex.vtk",tk);CHKERRQ(ierr);
224c4762a1bSJed Brown         ierr = PetscViewerFileSetName(viewer,fname);CHKERRQ(ierr);
225c4762a1bSJed Brown         ierr = DMView(celldm,viewer);CHKERRQ(ierr);
226c4762a1bSJed Brown         ierr = VecView(pfields[0],viewer);CHKERRQ(ierr);
227c4762a1bSJed Brown       }
228c4762a1bSJed Brown       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
229c4762a1bSJed Brown     }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   }
232c4762a1bSJed Brown   ierr = VecDestroy(&pfields[0]);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = PetscFree(pfields);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = DMDestroy(&celldm);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = DMDestroy(&swarm);CHKERRQ(ierr);
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown int main(int argc,char **args)
241c4762a1bSJed Brown {
242c4762a1bSJed Brown   PetscErrorCode ierr;
243c4762a1bSJed Brown   PetscInt ppcell = 1;
244c4762a1bSJed Brown   PetscInt meshtype = 0;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
247c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-ppcell",&ppcell,NULL);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-meshtype",&meshtype,NULL);CHKERRQ(ierr);
249*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(meshtype > 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"-meshtype <value> must be 0 or 1");
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   ierr = pic_advect(ppcell,meshtype);CHKERRQ(ierr);
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   ierr = PetscFinalize();
254c4762a1bSJed Brown   return ierr;
255c4762a1bSJed Brown }
256