xref: /petsc/src/dm/tutorials/ex21.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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   const PetscInt dim = 2;
17c4762a1bSJed Brown   DM celldm,swarm;
18c4762a1bSJed Brown   PetscInt tk,nt = 200;
19c4762a1bSJed Brown   PetscBool view = PETSC_FALSE;
20c4762a1bSJed Brown   Vec *pfields;
21c4762a1bSJed Brown   PetscReal minradius;
22c4762a1bSJed Brown   PetscReal dt;
23c4762a1bSJed Brown   PetscReal vel[] = { 1.0, 0.16 };
24c4762a1bSJed Brown   const char *fieldnames[] = { "phi" };
25c4762a1bSJed Brown   PetscViewer viewer;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   PetscFunctionBegin;
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nt",&nt,NULL));
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   /* Create the background cell DM */
32c4762a1bSJed Brown   if (meshtype == 0) { /* DA */
33c4762a1bSJed Brown     PetscInt nxy;
34c4762a1bSJed Brown     PetscInt dof = 1;
35c4762a1bSJed Brown     PetscInt stencil_width = 1;
36c4762a1bSJed Brown 
379566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMDA\n"));
38c4762a1bSJed Brown     nxy = 33;
399566063dSJacob 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));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch     PetscCall(DMDASetElementType(celldm,DMDA_ELEMENT_Q1));
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(celldm));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch     PetscCall(DMSetUp(celldm));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch     PetscCall(DMDASetUniformCoordinates(celldm,0.0,1.0,0.0,1.0,0.0,1.5));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown     minradius = 1.0/((PetscReal)(nxy-1));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"DA(minradius) %1.4e\n",(double)minradius));
52c4762a1bSJed Brown   }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   if (meshtype == 1){ /* PLEX */
55c4762a1bSJed Brown     DM distributedMesh = NULL;
56c4762a1bSJed Brown     PetscInt numComp[] = {1};
57c4762a1bSJed Brown     PetscInt numDof[] = {1,0,0}; /* vert, edge, cell */
58c4762a1bSJed Brown     PetscInt faces[]  = {1,1,1};
59c4762a1bSJed Brown     PetscInt numBC = 0;
60c4762a1bSJed Brown     PetscSection section;
61c4762a1bSJed Brown     Vec cellgeom = NULL;
62c4762a1bSJed Brown     Vec facegeom = NULL;
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mesh type: DMPLEX\n"));
659566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, &celldm));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown     /* Distribute mesh over processes */
689566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(celldm,0,NULL,&distributedMesh));
69c4762a1bSJed Brown     if (distributedMesh) {
709566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&celldm));
71c4762a1bSJed Brown       celldm = distributedMesh;
72c4762a1bSJed Brown     }
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(celldm));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(celldm,NULL,numComp,numDof,numBC,NULL,NULL,NULL,NULL,&section));
779566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(celldm,section));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch     PetscCall(DMSetUp(celldm));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */
829566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeGeometryFVM(celldm,&cellgeom,&facegeom));
839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetMinRadius(celldm,&minradius));
849566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"PLEX(minradius) %1.4e\n",(double)minradius));
859566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&cellgeom));
869566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&facegeom));
879566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* Create the DMSwarm */
919566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD,&swarm));
929566063dSJacob Faibussowitsch   PetscCall(DMSetType(swarm,DMSWARM));
939566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(swarm,dim));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Configure swarm to be of type PIC */
969566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(swarm,DMSWARM_PIC));
979566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(swarm,celldm));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* Register two scalar fields within the DMSwarm */
1009566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"phi",1,PETSC_REAL));
1019566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(swarm,"region",1,PETSC_REAL));
1029566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(swarm));
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /* Set initial local sizes of the DMSwarm with a buffer length of zero */
1059566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(swarm,4,0));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Insert swarm coordinates cell-wise */
1089566063dSJacob Faibussowitsch   /*PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell));*/
1099566063dSJacob Faibussowitsch   PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_SUBDIVISION,ppcell));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* Define initial conditions for th swarm fields "phi" and "region" */
112c4762a1bSJed Brown   {
113c4762a1bSJed Brown     PetscReal *s_coor,*s_phi,*s_region;
114c4762a1bSJed Brown     PetscInt npoints,p;
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(swarm,&npoints));
1179566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
1189566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,"phi",NULL,NULL,(void**)&s_phi));
1199566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,"region",NULL,NULL,(void**)&s_region));
120c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
121c4762a1bSJed Brown       PetscReal pos[2];
122c4762a1bSJed Brown       pos[0] = s_coor[2*p+0];
123c4762a1bSJed Brown       pos[1] = s_coor[2*p+1];
124c4762a1bSJed Brown 
125c4762a1bSJed Brown       s_region[p] = 1.0;
126c4762a1bSJed 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)));
127c4762a1bSJed Brown     }
1289566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,"region",NULL,NULL,(void**)&s_region));
1299566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,"phi",NULL,NULL,(void**)&s_phi));
1309566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Project initial value of phi onto the mesh */
1349566063dSJacob Faibussowitsch   PetscCall(DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_FALSE));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   if (view) {
137c4762a1bSJed Brown     /* View swarm all swarm fields using data type PETSC_REAL */
1389566063dSJacob Faibussowitsch     PetscCall(DMSwarmViewXDMF(swarm,"ic_dms.xmf"));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     /* View projected swarm field "phi" */
1419566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
1429566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
1439566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
144c4762a1bSJed Brown     if (meshtype == 0) { /* DA */
1459566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetName(viewer,"ic_dmda.vts"));
1469566063dSJacob Faibussowitsch       PetscCall(VecView(pfields[0],viewer));
147c4762a1bSJed Brown     }
148c4762a1bSJed Brown     if (meshtype == 1) { /* PLEX */
1499566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetName(viewer,"ic_dmplex.vtk"));
1509566063dSJacob Faibussowitsch       PetscCall(DMView(celldm,viewer));
1519566063dSJacob Faibussowitsch       PetscCall(VecView(pfields[0],viewer));
152c4762a1bSJed Brown     }
1539566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
154c4762a1bSJed Brown   }
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch   PetscCall(DMView(celldm,PETSC_VIEWER_STDOUT_WORLD));
1579566063dSJacob Faibussowitsch   PetscCall(DMView(swarm,PETSC_VIEWER_STDOUT_WORLD));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   dt = 0.5 * minradius / PetscSqrtReal(vel[0]*vel[0] + vel[1]*vel[1]);
160c4762a1bSJed Brown   for (tk=1; tk<=nt; tk++) {
161c4762a1bSJed Brown     PetscReal *s_coor;
162c4762a1bSJed Brown     PetscInt npoints,p;
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[step %D]\n",tk));
165c4762a1bSJed Brown     /* advect with analytic prescribed (constant) velocity field */
1669566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(swarm,&npoints));
1679566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
168c4762a1bSJed Brown     for (p=0; p<npoints; p++) {
169c4762a1bSJed Brown       s_coor[2*p+0] += dt * vel[0];
170c4762a1bSJed Brown       s_coor[2*p+1] += dt * vel[1];
171c4762a1bSJed Brown     }
1729566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(swarm,DMSwarmPICField_coor,NULL,NULL,(void**)&s_coor));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(swarm,PETSC_TRUE));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown     /* Ad-hoc cell filling algorithm */
177c4762a1bSJed Brown     /*
178c4762a1bSJed Brown        The injection frequency is chosen for default DA case.
179c4762a1bSJed Brown        They will likely not be appropriate for the general case using an unstructure PLEX mesh.
180c4762a1bSJed Brown     */
181c4762a1bSJed Brown     if (tk%10 == 0) {
182c4762a1bSJed Brown       PetscReal dx = 1.0/32.0;
183c4762a1bSJed Brown       PetscInt npoints_dir_x[] = { 32, 1 };
184c4762a1bSJed Brown       PetscReal min[2],max[2];
185c4762a1bSJed Brown 
186c4762a1bSJed Brown       min[0] = 0.5 * dx;  max[0] = 0.5 * dx + 31.0 * dx;
187c4762a1bSJed Brown       min[1] = 0.5 * dx;  max[1] = 0.5 * dx;
1889566063dSJacob Faibussowitsch       PetscCall(DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_x,ADD_VALUES));
189c4762a1bSJed Brown     }
190c4762a1bSJed Brown     if (tk%2 == 0) {
191c4762a1bSJed Brown       PetscReal dx = 1.0/32.0;
192c4762a1bSJed Brown       PetscInt npoints_dir_y[] = { 2, 31 };
193c4762a1bSJed Brown       PetscReal min[2],max[2];
194c4762a1bSJed Brown 
195c4762a1bSJed Brown       min[0] = 0.05 * dx; max[0] = 0.5 * dx;
196c4762a1bSJed Brown       min[1] = 0.5 * dx;  max[1] = 0.5 * dx + 31.0 * dx;
1979566063dSJacob Faibussowitsch       PetscCall(DMSwarmSetPointsUniformCoordinates(swarm,min,max,npoints_dir_y,ADD_VALUES));
198c4762a1bSJed Brown     }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown     /* Project swarm field "phi" onto the cell DM */
2019566063dSJacob Faibussowitsch     PetscCall(DMSwarmProjectFields(swarm,1,fieldnames,&pfields,PETSC_TRUE));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown     if (view) {
204c4762a1bSJed Brown       PetscViewer viewer;
205c4762a1bSJed Brown       char fname[PETSC_MAX_PATH_LEN];
206c4762a1bSJed Brown 
207c4762a1bSJed Brown       /* View swarm fields */
2089566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dms.xmf",tk));
2099566063dSJacob Faibussowitsch       PetscCall(DMSwarmViewXDMF(swarm,fname));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown       /* View projected field */
2129566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
2139566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer,PETSCVIEWERVTK));
2149566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_WRITE));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown       if (meshtype == 0) { /* DA */
2179566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmda.vts",tk));
2189566063dSJacob Faibussowitsch         PetscCall(PetscViewerFileSetName(viewer,fname));
2199566063dSJacob Faibussowitsch         PetscCall(VecView(pfields[0],viewer));
220c4762a1bSJed Brown       }
221c4762a1bSJed Brown       if (meshtype == 1) { /* PLEX */
2229566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fname,PETSC_MAX_PATH_LEN-1,"step%.4D_dmplex.vtk",tk));
2239566063dSJacob Faibussowitsch         PetscCall(PetscViewerFileSetName(viewer,fname));
2249566063dSJacob Faibussowitsch         PetscCall(DMView(celldm,viewer));
2259566063dSJacob Faibussowitsch         PetscCall(VecView(pfields[0],viewer));
226c4762a1bSJed Brown       }
2279566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
228c4762a1bSJed Brown     }
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   }
2319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pfields[0]));
2329566063dSJacob Faibussowitsch   PetscCall(PetscFree(pfields));
2339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&celldm));
2349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&swarm));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown int main(int argc,char **args)
240c4762a1bSJed Brown {
241c4762a1bSJed Brown   PetscInt ppcell = 1;
242c4762a1bSJed Brown   PetscInt meshtype = 0;
243c4762a1bSJed Brown 
2449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
2459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-ppcell",&ppcell,NULL));
2469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-meshtype",&meshtype,NULL));
247*08401ef6SPierre Jolivet   PetscCheck(meshtype <= 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"-meshtype <value> must be 0 or 1");
248c4762a1bSJed Brown 
2499566063dSJacob Faibussowitsch   PetscCall(pic_advect(ppcell,meshtype));
250c4762a1bSJed Brown 
2519566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
252b122ec5aSJacob Faibussowitsch   return 0;
253c4762a1bSJed Brown }
254