xref: /petsc/src/dm/tutorials/ex21.c (revision 9566063d113dddea24716c546802770db7481bc0)
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;
29*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL));
30*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nt",&nt,NULL));
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 
38*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMDA\n"));
39c4762a1bSJed Brown     nxy = 33;
40*9566063dSJacob Faibussowitsch     PetscCall(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));
41c4762a1bSJed Brown 
42*9566063dSJacob Faibussowitsch     PetscCall(DMDASetElementType(celldm,DMDA_ELEMENT_Q1));
43c4762a1bSJed Brown 
44*9566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(celldm));
45c4762a1bSJed Brown 
46*9566063dSJacob Faibussowitsch     PetscCall(DMSetUp(celldm));
47c4762a1bSJed Brown 
48*9566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(celldm,0.0,1.0,0.0,1.0,0.0,1.5));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown     minradius = 1.0/((PetscReal)(nxy-1));
51c4762a1bSJed Brown 
52*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DA(minradius) %1.4e\n",(double)minradius));
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 
65*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMPLEX\n"));
66*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, &celldm));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown     /* Distribute mesh over processes */
69*9566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh));
70c4762a1bSJed Brown     if (distributedMesh) {
71*9566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&celldm));
72c4762a1bSJed Brown       celldm = distributedMesh;
73c4762a1bSJed Brown     }
74c4762a1bSJed Brown 
75*9566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(celldm));
76c4762a1bSJed Brown 
77*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section));
78*9566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(celldm,section));
79c4762a1bSJed Brown 
80*9566063dSJacob Faibussowitsch     PetscCall(DMSetUp(celldm));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown     /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */
83*9566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeGeometryFVM(celldm,&cellgeom,&facegeom));
84*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetMinRadius(celldm,&minradius));
85*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"PLEX(minradius) %1.4e\n",(double)minradius));
86*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cellgeom));
87*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&facegeom));
88*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Create the DMSwarm */
92*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm));
93*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm,DMSWARM));
94*9566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm,dim));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
97*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC));
98*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm,celldm));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
101*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"phi",1,PETSC_REAL));
102*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"region",1,PETSC_REAL));
103*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
106*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm,4,0));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
109*9566063dSJacob Faibussowitsch   /*PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell));*/
110*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,ppcell));
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 
117*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(swarm,&npoints));
118*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
119*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,"phi",NULL,NULL,(void**)&s_phi));
120*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,"region",NULL,NULL,(void**)&s_region));
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     }
129*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,"region",NULL,NULL,(void**)&s_region));
130*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,"phi",NULL,NULL,(void**)&s_phi));
131*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   /* Project initial value of phi onto the mesh */
135*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_FALSE));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   if (view) {
138c4762a1bSJed Brown     /* View swarm all swarm fields using data type PETSC_REAL */
139*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmViewXDMF(swarm,"ic_dms.xmf"));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown     /* View projected swarm field "phi" */
142*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
143*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
144*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
145c4762a1bSJed Brown     if (meshtype == 0) { /* DA */
146*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetName(viewer,"ic_dmda.vts"));
147*9566063dSJacob Faibussowitsch       PetscCall(VecView(pfields[0],viewer));
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown     if (meshtype == 1) { /* PLEX */
150*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetName(viewer,"ic_dmplex.vtk"));
151*9566063dSJacob Faibussowitsch       PetscCall(DMView(celldm,viewer));
152*9566063dSJacob Faibussowitsch       PetscCall(VecView(pfields[0],viewer));
153c4762a1bSJed Brown     }
154*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown 
157*9566063dSJacob Faibussowitsch   PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
158*9566063dSJacob Faibussowitsch   PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
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 
165*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[step %D]\n",tk));
166c4762a1bSJed Brown     /* advect with analytic prescribed (constant) velocity field */
167*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(swarm,&npoints));
168*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
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     }
173*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
174c4762a1bSJed Brown 
175*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(swarm,PETSC_TRUE));
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;
189*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_x,ADD_VALUES));
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;
198*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_y,ADD_VALUES));
199c4762a1bSJed Brown     }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown     /* Project swarm field "phi" onto the cell DM */
202*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_TRUE));
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 */
209*9566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dms.xmf",tk));
210*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmViewXDMF(swarm,fname));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown       /* View projected field */
213*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
214*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
215*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown       if (meshtype == 0) { /* DA */
218*9566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmda.vts",tk));
219*9566063dSJacob Faibussowitsch         PetscCall(PetscViewerFileSetName(viewer,fname));
220*9566063dSJacob Faibussowitsch         PetscCall(VecView(pfields[0],viewer));
221c4762a1bSJed Brown       }
222c4762a1bSJed Brown       if (meshtype == 1) { /* PLEX */
223*9566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmplex.vtk",tk));
224*9566063dSJacob Faibussowitsch         PetscCall(PetscViewerFileSetName(viewer,fname));
225*9566063dSJacob Faibussowitsch         PetscCall(DMView(celldm,viewer));
226*9566063dSJacob Faibussowitsch         PetscCall(VecView(pfields[0],viewer));
227c4762a1bSJed Brown       }
228*9566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
229c4762a1bSJed Brown     }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   }
232*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pfields[0]));
233*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(pfields));
234*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
235*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
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 
246*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
247*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ppcell",&ppcell,NULL));
248*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-meshtype",&meshtype,NULL));
2492c71b3e2SJacob Faibussowitsch   PetscCheckFalse(meshtype > 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"-meshtype <value> must be 0 or 1");
250c4762a1bSJed Brown 
251*9566063dSJacob Faibussowitsch   PetscCall(pic_advect(ppcell,meshtype));
252c4762a1bSJed Brown 
253*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
254b122ec5aSJacob Faibussowitsch   return 0;
255c4762a1bSJed Brown }
256