1145b44c9SPierre Jolivet static char help[] = "DMSwarm-PIC demonstrator of advecting points within cell DM defined by a DA or PLEX object \n\
2c4762a1bSJed Brown Options: \n\
3c4762a1bSJed Brown -ppcell : Number of times to sub-divide the reference cell when layout the initial particle coordinates \n\
4c4762a1bSJed Brown -meshtype : 0 ==> DA , 1 ==> PLEX \n\
5c4762a1bSJed Brown -nt : Number of timestep to perform \n\
6c4762a1bSJed Brown -view : Write out initial condition and time dependent data \n";
7c4762a1bSJed Brown
8c4762a1bSJed Brown #include <petsc.h>
9c4762a1bSJed Brown #include <petscdm.h>
10c4762a1bSJed Brown #include <petscdmda.h>
11c4762a1bSJed Brown #include <petscdmswarm.h>
12c4762a1bSJed Brown
pic_advect(PetscInt ppcell,PetscInt meshtype)13d71ae5a4SJacob Faibussowitsch PetscErrorCode pic_advect(PetscInt ppcell, PetscInt meshtype)
14d71ae5a4SJacob Faibussowitsch {
15c4762a1bSJed Brown const PetscInt dim = 2;
16c4762a1bSJed Brown DM celldm, swarm;
17c4762a1bSJed Brown PetscInt tk, nt = 200;
18c4762a1bSJed Brown PetscBool view = PETSC_FALSE;
19fc3eb83cSMatthew G. Knepley Vec pfields[1];
20c4762a1bSJed Brown PetscReal minradius;
21c4762a1bSJed Brown PetscReal dt;
22c4762a1bSJed Brown PetscReal vel[] = {1.0, 0.16};
23c4762a1bSJed Brown const char *fieldnames[] = {"phi"};
24c4762a1bSJed Brown PetscViewer viewer;
25c4762a1bSJed Brown
26c4762a1bSJed Brown PetscFunctionBegin;
279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view", &view, NULL));
289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nt", &nt, NULL));
29c4762a1bSJed Brown
30c4762a1bSJed Brown /* Create the background cell DM */
31c4762a1bSJed Brown if (meshtype == 0) { /* DA */
32c4762a1bSJed Brown PetscInt nxy;
33c4762a1bSJed Brown PetscInt dof = 1;
34c4762a1bSJed Brown PetscInt stencil_width = 1;
35c4762a1bSJed Brown
369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Mesh type: DMDA\n"));
37c4762a1bSJed Brown nxy = 33;
389566063dSJacob 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));
39c4762a1bSJed Brown
409566063dSJacob Faibussowitsch PetscCall(DMDASetElementType(celldm, DMDA_ELEMENT_Q1));
41c4762a1bSJed Brown
429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm));
43c4762a1bSJed Brown
449566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm));
45c4762a1bSJed Brown
469566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(celldm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.5));
47c4762a1bSJed Brown
48c4762a1bSJed Brown minradius = 1.0 / ((PetscReal)(nxy - 1));
49c4762a1bSJed Brown
509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DA(minradius) %1.4e\n", (double)minradius));
51c4762a1bSJed Brown }
52c4762a1bSJed Brown
53c4762a1bSJed Brown if (meshtype == 1) { /* PLEX */
54c4762a1bSJed Brown DM distributedMesh = NULL;
55c4762a1bSJed Brown PetscInt numComp[] = {1};
56c4762a1bSJed Brown PetscInt numDof[] = {1, 0, 0}; /* vert, edge, cell */
57c4762a1bSJed Brown PetscInt faces[] = {1, 1, 1};
58c4762a1bSJed Brown PetscInt numBC = 0;
59c4762a1bSJed Brown PetscSection section;
60c4762a1bSJed Brown Vec cellgeom = NULL;
61c4762a1bSJed Brown Vec facegeom = NULL;
62c4762a1bSJed Brown
639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Mesh type: DMPLEX\n"));
6442108689Sksagiyam PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_TRUE, faces, NULL, NULL, PETSC_TRUE, 0, PETSC_TRUE, &celldm));
65c4762a1bSJed Brown
66c4762a1bSJed Brown /* Distribute mesh over processes */
679566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(celldm, 0, NULL, &distributedMesh));
68c4762a1bSJed Brown if (distributedMesh) {
699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm));
70c4762a1bSJed Brown celldm = distributedMesh;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown
739566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(celldm));
74c4762a1bSJed Brown
759566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(celldm, NULL, numComp, numDof, numBC, NULL, NULL, NULL, NULL, §ion));
769566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(celldm, section));
77c4762a1bSJed Brown
789566063dSJacob Faibussowitsch PetscCall(DMSetUp(celldm));
79c4762a1bSJed Brown
80c4762a1bSJed Brown /* Calling DMPlexComputeGeometryFVM() generates the value returned by DMPlexGetMinRadius() */
819566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGeometryFVM(celldm, &cellgeom, &facegeom));
829566063dSJacob Faibussowitsch PetscCall(DMPlexGetMinRadius(celldm, &minradius));
839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "PLEX(minradius) %1.4e\n", (double)minradius));
849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cellgeom));
859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&facegeom));
869566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion));
87c4762a1bSJed Brown }
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* Create the DMSwarm */
909566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &swarm));
919566063dSJacob Faibussowitsch PetscCall(DMSetType(swarm, DMSWARM));
929566063dSJacob Faibussowitsch PetscCall(DMSetDimension(swarm, dim));
93c4762a1bSJed Brown
94c4762a1bSJed Brown /* Configure swarm to be of type PIC */
959566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(swarm, DMSWARM_PIC));
969566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(swarm, celldm));
97c4762a1bSJed Brown
98c4762a1bSJed Brown /* Register two scalar fields within the DMSwarm */
999566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "phi", 1, PETSC_REAL));
1009566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(swarm, "region", 1, PETSC_REAL));
1019566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(swarm));
102c4762a1bSJed Brown
103c4762a1bSJed Brown /* Set initial local sizes of the DMSwarm with a buffer length of zero */
1049566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(swarm, 4, 0));
105c4762a1bSJed Brown
106c4762a1bSJed Brown /* Insert swarm coordinates cell-wise */
1079566063dSJacob Faibussowitsch /*PetscCall(DMSwarmInsertPointsUsingCellDM(swarm,DMSWARMPIC_LAYOUT_REGULAR,ppcell));*/
1089566063dSJacob Faibussowitsch PetscCall(DMSwarmInsertPointsUsingCellDM(swarm, DMSWARMPIC_LAYOUT_SUBDIVISION, ppcell));
109c4762a1bSJed Brown
110c4762a1bSJed Brown /* Define initial conditions for th swarm fields "phi" and "region" */
111c4762a1bSJed Brown {
112c4762a1bSJed Brown PetscReal *s_coor, *s_phi, *s_region;
113c4762a1bSJed Brown PetscInt npoints, p;
114c4762a1bSJed Brown
1159566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(swarm, &npoints));
1169566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&s_coor));
1179566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, "phi", NULL, NULL, (void **)&s_phi));
1189566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, "region", NULL, NULL, (void **)&s_region));
119c4762a1bSJed Brown for (p = 0; p < npoints; p++) {
120c4762a1bSJed Brown PetscReal pos[2];
121c4762a1bSJed Brown pos[0] = s_coor[2 * p + 0];
122c4762a1bSJed Brown pos[1] = s_coor[2 * p + 1];
123c4762a1bSJed Brown
124c4762a1bSJed Brown s_region[p] = 1.0;
125c4762a1bSJed 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)));
126c4762a1bSJed Brown }
1279566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, "region", NULL, NULL, (void **)&s_region));
1289566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, "phi", NULL, NULL, (void **)&s_phi));
1299566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&s_coor));
130c4762a1bSJed Brown }
131c4762a1bSJed Brown
132c4762a1bSJed Brown /* Project initial value of phi onto the mesh */
133fc3eb83cSMatthew G. Knepley PetscCall(DMCreateGlobalVector(celldm, &pfields[0]));
134953494dbSMatthew G. Knepley PetscCall(DMSwarmProjectFields(swarm, NULL, 1, fieldnames, pfields, SCATTER_FORWARD));
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
16463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[step %" PetscInt_FMT "]\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
1869371c9d4SSatish Balay min[0] = 0.5 * dx;
1879371c9d4SSatish Balay max[0] = 0.5 * dx + 31.0 * dx;
1889371c9d4SSatish Balay min[1] = 0.5 * dx;
1899371c9d4SSatish Balay max[1] = 0.5 * dx;
1909566063dSJacob Faibussowitsch PetscCall(DMSwarmSetPointsUniformCoordinates(swarm, min, max, npoints_dir_x, ADD_VALUES));
191c4762a1bSJed Brown }
192c4762a1bSJed Brown if (tk % 2 == 0) {
193c4762a1bSJed Brown PetscReal dx = 1.0 / 32.0;
194c4762a1bSJed Brown PetscInt npoints_dir_y[] = {2, 31};
195c4762a1bSJed Brown PetscReal min[2], max[2];
196c4762a1bSJed Brown
1979371c9d4SSatish Balay min[0] = 0.05 * dx;
1989371c9d4SSatish Balay max[0] = 0.5 * dx;
1999371c9d4SSatish Balay min[1] = 0.5 * dx;
2009371c9d4SSatish Balay max[1] = 0.5 * dx + 31.0 * dx;
2019566063dSJacob Faibussowitsch PetscCall(DMSwarmSetPointsUniformCoordinates(swarm, min, max, npoints_dir_y, ADD_VALUES));
202c4762a1bSJed Brown }
203c4762a1bSJed Brown
204c4762a1bSJed Brown /* Project swarm field "phi" onto the cell DM */
205953494dbSMatthew G. Knepley PetscCall(DMSwarmProjectFields(swarm, NULL, 1, fieldnames, pfields, SCATTER_FORWARD));
206c4762a1bSJed Brown
207c4762a1bSJed Brown if (view) {
208c4762a1bSJed Brown PetscViewer viewer;
209c4762a1bSJed Brown char fname[PETSC_MAX_PATH_LEN];
210c4762a1bSJed Brown
211c4762a1bSJed Brown /* View swarm fields */
21263a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, PETSC_MAX_PATH_LEN - 1, "step%.4" PetscInt_FMT "_dms.xmf", tk));
2139566063dSJacob Faibussowitsch PetscCall(DMSwarmViewXDMF(swarm, fname));
214c4762a1bSJed Brown
215c4762a1bSJed Brown /* View projected field */
2169566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer));
2179566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERVTK));
2189566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
219c4762a1bSJed Brown
220c4762a1bSJed Brown if (meshtype == 0) { /* DA */
22163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, PETSC_MAX_PATH_LEN - 1, "step%.4" PetscInt_FMT "_dmda.vts", tk));
2229566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, fname));
2239566063dSJacob Faibussowitsch PetscCall(VecView(pfields[0], viewer));
224c4762a1bSJed Brown }
225c4762a1bSJed Brown if (meshtype == 1) { /* PLEX */
22663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fname, PETSC_MAX_PATH_LEN - 1, "step%.4" PetscInt_FMT "_dmplex.vtk", tk));
2279566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, fname));
2289566063dSJacob Faibussowitsch PetscCall(DMView(celldm, viewer));
2299566063dSJacob Faibussowitsch PetscCall(VecView(pfields[0], viewer));
230c4762a1bSJed Brown }
2319566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
232c4762a1bSJed Brown }
233c4762a1bSJed Brown }
2349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pfields[0]));
2359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&celldm));
2369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&swarm));
2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown
main(int argc,char ** args)240d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
241d71ae5a4SJacob Faibussowitsch {
242c4762a1bSJed Brown PetscInt ppcell = 1;
243c4762a1bSJed Brown PetscInt meshtype = 0;
244c4762a1bSJed Brown
245327415f7SBarry Smith PetscFunctionBeginUser;
246*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
2479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ppcell", &ppcell, NULL));
2489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-meshtype", &meshtype, NULL));
24908401ef6SPierre Jolivet PetscCheck(meshtype <= 1, PETSC_COMM_WORLD, PETSC_ERR_USER, "-meshtype <value> must be 0 or 1");
250c4762a1bSJed Brown
2519566063dSJacob Faibussowitsch PetscCall(pic_advect(ppcell, meshtype));
252c4762a1bSJed Brown
2539566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
254b122ec5aSJacob Faibussowitsch return 0;
255c4762a1bSJed Brown }
256