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