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