xref: /petsc/src/dm/impls/swarm/tests/ex7.c (revision 1511cd715a1f0c8d257549c5ebe5cee9c6feed4d)
1 static char help[] = "Example program demonstrating projection between particle and finite element spaces using OpenMP in 2D cylindrical coordinates\n";
2 
3 #include "petscdmplex.h"
4 #include "petscds.h"
5 #include "petscdmswarm.h"
6 #include "petscksp.h"
7 #include <petsc/private/petscimpl.h>
8 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
9 #include <omp.h>
10 #endif
11 
12 typedef struct {
13   Mat MpTrans;
14   Mat Mp;
15   Vec ff;
16   Vec uu;
17 } MatShellCtx;
18 
19 PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy) {
20   MatShellCtx *matshellctx;
21 
22   PetscFunctionBeginUser;
23   PetscCall(MatShellGetContext(MtM, &matshellctx));
24   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
25   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
26   PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
27   PetscFunctionReturn(0);
28 }
29 
30 PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz) {
31   MatShellCtx *matshellctx;
32 
33   PetscFunctionBeginUser;
34   PetscCall(MatShellGetContext(MtM, &matshellctx));
35   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
36   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
37   PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
38   PetscFunctionReturn(0);
39 }
40 
41 PetscErrorCode createSwarm(const DM dm, DM *sw) {
42   PetscInt Nc = 1, dim = 2;
43 
44   PetscFunctionBeginUser;
45   PetscCall(DMCreate(PETSC_COMM_SELF, sw));
46   PetscCall(DMSetType(*sw, DMSWARM));
47   PetscCall(DMSetDimension(*sw, dim));
48   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
49   PetscCall(DMSwarmSetCellDM(*sw, dm));
50   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR));
51   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
52   PetscCall(DMSetFromOptions(*sw));
53   PetscFunctionReturn(0);
54 }
55 
56 PetscErrorCode gridToParticles(const DM dm, DM sw, PetscReal *moments, Vec rhs, Mat M_p) {
57   PetscBool     is_lsqr;
58   KSP           ksp;
59   Mat           PM_p = NULL, MtM, D;
60   Vec           ff;
61   PetscInt      Np, timestep = 0, bs, N, M, nzl;
62   PetscReal     time = 0.0;
63   PetscDataType dtype;
64   MatShellCtx  *matshellctx;
65 
66   PetscFunctionBeginUser;
67   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
68   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
69   PetscCall(KSPSetFromOptions(ksp));
70   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
71   if (!is_lsqr) {
72     PetscCall(MatGetLocalSize(M_p, &M, &N));
73     if (N > M) {
74       PC pc;
75       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") -- skip revert to lsqr\n", M, N));
76       is_lsqr = PETSC_TRUE;
77       PetscCall(KSPSetType(ksp, KSPLSQR));
78       PetscCall(KSPGetPC(ksp, &pc));
79       PetscCall(PCSetType(pc, PCNONE)); // could put in better solver -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero
80     } else {
81       PetscCall(PetscNew(&matshellctx));
82       PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
83       PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
84       matshellctx->Mp = M_p;
85       PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
86       PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
87       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
88       PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
89       for (int i = 0; i < N; i++) {
90         const PetscScalar *vals;
91         const PetscInt    *cols;
92         PetscScalar        dot = 0;
93         PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
94         for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
95         PetscCheck(dot != 0.0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Row %d is empty", i);
96         PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
97       }
98       PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
99       PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
100       PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl);
101       PetscCall(KSPSetOperators(ksp, MtM, D));
102       PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
103       PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
104       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
105     }
106   }
107   if (is_lsqr) {
108     PC        pc;
109     PetscBool is_bjac;
110     PetscCall(KSPGetPC(ksp, &pc));
111     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &is_bjac));
112     if (is_bjac) {
113       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
114       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
115     } else {
116       PetscCall(KSPSetOperators(ksp, M_p, M_p));
117     }
118   }
119   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!!
120   if (!is_lsqr) {
121     PetscCall(KSPSolve(ksp, rhs, matshellctx->uu));
122     PetscCall(MatMult(M_p, matshellctx->uu, ff));
123     PetscCall(MatDestroy(&matshellctx->MpTrans));
124     PetscCall(VecDestroy(&matshellctx->ff));
125     PetscCall(VecDestroy(&matshellctx->uu));
126     PetscCall(MatDestroy(&D));
127     PetscCall(MatDestroy(&MtM));
128     PetscCall(PetscFree(matshellctx));
129   } else {
130     PetscCall(KSPSolveTranspose(ksp, rhs, ff));
131   }
132   PetscCall(KSPDestroy(&ksp));
133   /* Visualize particle field */
134   PetscCall(DMSetOutputSequenceNumber(sw, timestep, time));
135   PetscCall(VecViewFromOptions(ff, NULL, "-weights_view"));
136   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
137 
138   /* compute energy */
139   if (moments) {
140     PetscReal *wq, *coords;
141     PetscCall(DMSwarmGetLocalSize(sw, &Np));
142     PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
143     PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
144     moments[0] = moments[1] = moments[2] = 0;
145     for (int p = 0; p < Np; p++) {
146       moments[0] += wq[p];
147       moments[1] += wq[p] * coords[p * 2 + 0]; // x-momentum
148       moments[2] += wq[p] * (PetscSqr(coords[p * 2 + 0]) + PetscSqr(coords[p * 2 + 1]));
149     }
150     PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
151     PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
152   }
153   PetscCall(MatDestroy(&PM_p));
154   PetscFunctionReturn(0);
155 }
156 
157 PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt Np, const PetscInt a_tid, const PetscInt dim, const PetscInt target, const PetscReal xx[], const PetscReal yy[], const PetscReal a_wp[], Vec rho, Mat *Mp_out) {
158   PetscBool     removePoints = PETSC_TRUE;
159   PetscReal    *wq, *coords;
160   PetscDataType dtype;
161   Mat           M_p;
162   Vec           ff;
163   PetscInt      bs, p, zero = 0;
164 
165   PetscFunctionBeginUser;
166   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
167   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
168   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
169   for (p = 0; p < Np; p++) {
170     coords[p * 2 + 0] = xx[p];
171     coords[p * 2 + 1] = yy[p];
172     wq[p]             = a_wp[p];
173   }
174   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
175   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
176   PetscCall(DMSwarmMigrate(sw, removePoints));
177   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
178   if (a_tid == target) PetscCall(DMViewFromOptions(sw, NULL, "-swarm_view"));
179 
180   /* Project particles to field */
181   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
182   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
183   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
184 
185   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!!
186   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
187   PetscCall(MatMultTranspose(M_p, ff, rho));
188   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
189 
190   /* Visualize mesh field */
191   if (a_tid == target) PetscCall(VecViewFromOptions(rho, NULL, "-rho_view"));
192   // output
193   *Mp_out = M_p;
194 
195   PetscFunctionReturn(0);
196 }
197 static PetscErrorCode maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscScalar *u) {
198   PetscInt  i;
199   PetscReal v2 = 0, theta = 2 * kt_m; /* theta = 2kT/mc^2 */
200 
201   PetscFunctionBegin;
202   /* compute the exponents, v^2 */
203   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
204   /* evaluate the Maxwellian */
205   u[0] = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)) * 2. * PETSC_PI * x[1]; // radial term for 2D axi-sym.
206   PetscFunctionReturn(0);
207 }
208 
209 #define MAX_NUM_THRDS 12
210 PetscErrorCode go() {
211   DM        dm_t[MAX_NUM_THRDS], sw_t[MAX_NUM_THRDS];
212   PetscFE   fe;
213   PetscInt  dim = 2, Nc = 1, i, faces[3];
214   PetscInt  Np[2] = {10, 10}, Np2[2], field = 0, target = 0, Np_t[MAX_NUM_THRDS];
215   PetscReal moments_0[3], moments_1[3], vol = 1;
216   PetscReal lo[3] = {-5, 0, -5}, hi[3] = {5, 5, 5}, h[3], hp[3], *xx_t[MAX_NUM_THRDS], *yy_t[MAX_NUM_THRDS], *wp_t[MAX_NUM_THRDS];
217   Vec       rho_t[MAX_NUM_THRDS], rhs_t[MAX_NUM_THRDS];
218   Mat       M_p_t[MAX_NUM_THRDS];
219 #if defined PETSC_USE_LOG
220   PetscLogStage stage;
221   PetscLogEvent swarm_create_ev, solve_ev, solve_loop_ev;
222 #endif
223 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
224   PetscInt numthreads = PetscNumOMPThreads;
225 #else
226   PetscInt numthreads = 1;
227 #endif
228 
229   PetscFunctionBeginUser;
230 #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
231   PetscCheck(numthreads <= MAX_NUM_THRDS, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %" PetscInt_FMT " > %d", numthreads, MAX_NUM_THRDS);
232   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %" PetscInt_FMT " > %d", numthreads, MAX_NUM_THRDS);
233 #endif
234   if (target >= numthreads) target = numthreads - 1;
235   PetscCall(PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev));
236   PetscCall(PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev));
237   PetscCall(PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev));
238   PetscCall(PetscLogStageRegister("Solve", &stage));
239   i = dim;
240   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL));
241   i = dim;
242   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-np", Np, &i, NULL));
243   /* Create thread meshes */
244   for (int tid = 0; tid < numthreads; tid++) {
245     // setup mesh dm_t, could use PETSc's Landau create velocity space mesh here to get dm_t[tid]
246     PetscCall(DMCreate(PETSC_COMM_SELF, &dm_t[tid]));
247     PetscCall(DMSetType(dm_t[tid], DMPLEX));
248     PetscCall(DMSetFromOptions(dm_t[tid]));
249     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe));
250     PetscCall(PetscFESetFromOptions(fe));
251     PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
252     PetscCall(DMSetField(dm_t[tid], field, NULL, (PetscObject)fe));
253     PetscCall(DMCreateDS(dm_t[tid]));
254     PetscCall(PetscFEDestroy(&fe));
255     // helper vectors
256     PetscCall(DMCreateGlobalVector(dm_t[tid], &rho_t[tid]));
257     PetscCall(DMCreateGlobalVector(dm_t[tid], &rhs_t[tid]));
258     // this mimics application code
259     PetscCall(DMGetBoundingBox(dm_t[tid], lo, hi));
260     if (tid == target) {
261       PetscCall(DMViewFromOptions(dm_t[tid], NULL, "-dm_view"));
262       for (i = 0, vol = 1; i < dim; i++) {
263         h[i]  = (hi[i] - lo[i]) / faces[i];
264         hp[i] = (hi[i] - lo[i]) / Np[i];
265         vol *= (hi[i] - lo[i]);
266         PetscCall(PetscInfo(dm_t[tid], " lo = %g hi = %g n = %" PetscInt_FMT " h = %g hp = %g\n", (double)lo[i], (double)hi[i], faces[i], (double)h[i], (double)hp[i]));
267       }
268     }
269   }
270   // prepare particle data for problems. This mimics application code
271   PetscCall(PetscLogEventBegin(swarm_create_ev, 0, 0, 0, 0));
272   Np2[0] = Np[0];
273   Np2[1] = Np[1];
274   for (int tid = 0; tid < numthreads; tid++) { // change size of particle list a little
275     Np_t[tid] = Np2[0] * Np2[1];
276     PetscCall(PetscMalloc3(Np_t[tid], &xx_t[tid], Np_t[tid], &yy_t[tid], Np_t[tid], &wp_t[tid]));
277     if (tid == target) moments_0[0] = moments_0[1] = moments_0[2] = 0;
278     for (int pi = 0, pp = 0; pi < Np2[0]; pi++) {
279       for (int pj = 0; pj < Np2[1]; pj++, pp++) {
280         xx_t[tid][pp] = lo[0] + hp[0] / 2. + pi * hp[0];
281         yy_t[tid][pp] = lo[1] + hp[1] / 2. + pj * hp[1];
282         {
283           PetscReal x[] = {xx_t[tid][pp], yy_t[tid][pp]};
284           PetscCall(maxwellian(2, x, 1.0, vol / (PetscReal)Np_t[tid], &wp_t[tid][pp]));
285         }
286         if (tid == target) { //energy_0 += wp_t[tid][pp]*(PetscSqr(xx_t[tid][pp])+PetscSqr(yy_t[tid][pp]));
287           moments_0[0] += wp_t[tid][pp];
288           moments_0[1] += wp_t[tid][pp] * xx_t[tid][pp]; // x-momentum
289           moments_0[2] += wp_t[tid][pp] * (PetscSqr(xx_t[tid][pp]) + PetscSqr(yy_t[tid][pp]));
290         }
291       }
292     }
293     Np2[0]++;
294     Np2[1]++;
295   }
296   PetscCall(PetscLogEventEnd(swarm_create_ev, 0, 0, 0, 0));
297   PetscCall(PetscLogEventBegin(solve_ev, 0, 0, 0, 0));
298   /* Create particle swarm */
299   PetscPragmaOMP(parallel for)
300   for (int tid=0; tid<numthreads; tid++) {
301     PetscCallAbort(PETSC_COMM_SELF, createSwarm(dm_t[tid], &sw_t[tid]));
302   }
303   PetscPragmaOMP(parallel for)
304   for (int tid=0; tid<numthreads; tid++) {
305     PetscCallAbort(PETSC_COMM_SELF, particlesToGrid(dm_t[tid], sw_t[tid], Np_t[tid], tid, dim, target, xx_t[tid], yy_t[tid], wp_t[tid], rho_t[tid], &M_p_t[tid]));
306   }
307   /* Project field to particles */
308   /*   This gives f_p = M_p^+ M f */
309   PetscPragmaOMP(parallel for)
310   for (int tid=0; tid<numthreads; tid++) {
311     PetscCallAbort(PETSC_COMM_SELF, VecCopy(rho_t[tid], rhs_t[tid])); /* Identity: M^1 M rho */
312   }
313   PetscPragmaOMP(parallel for)
314   for (int tid=0; tid<numthreads; tid++) {
315     PetscCallAbort(PETSC_COMM_SELF, gridToParticles(dm_t[tid], sw_t[tid], (tid == target) ? moments_1 : NULL, rhs_t[tid], M_p_t[tid]));
316   }
317   /* Cleanup */
318   for (int tid = 0; tid < numthreads; tid++) {
319     PetscCall(MatDestroy(&M_p_t[tid]));
320     PetscCall(DMDestroy(&sw_t[tid]));
321   }
322   PetscCall(PetscLogEventEnd(solve_ev, 0, 0, 0, 0));
323   //
324   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Total number density: %20.12e (%20.12e); x-momentum = %g (%g); energy = %g error = %e, %" PetscInt_FMT " particles. Use %" PetscInt_FMT " threads\n", (double)moments_1[0], (double)moments_0[0], (double)moments_1[1], (double)moments_0[1], (double)moments_1[2], (double)((moments_1[2] - moments_0[2]) / moments_0[2]), Np[0] * Np[1], numthreads));
325   /* Cleanup */
326   for (int tid = 0; tid < numthreads; tid++) {
327     PetscCall(VecDestroy(&rho_t[tid]));
328     PetscCall(VecDestroy(&rhs_t[tid]));
329     PetscCall(DMDestroy(&dm_t[tid]));
330     PetscCall(PetscFree3(xx_t[tid], yy_t[tid], wp_t[tid]));
331   }
332   PetscFunctionReturn(0);
333 }
334 
335 int main(int argc, char **argv) {
336   PetscFunctionBeginUser;
337   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
338   PetscCall(go());
339   PetscCall(PetscFinalize());
340   return 0;
341 }
342 
343 /*TEST
344 
345   build:
346     requires: !complex
347 
348   test:
349     suffix: 0
350     requires: double triangle
351     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 10 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -ftop_ksp_type lsqr -ftop_pc_type none -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
352     filter: grep -v DM_ | grep -v atomic
353 
354   test:
355     suffix: 1
356     requires: double triangle
357     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 10 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -dm_plex_hash_location -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
358     filter: grep -v DM_ | grep -v atomic
359 
360   test:
361     suffix: 2
362     requires: double triangle
363     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 10 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -dm_plex_hash_location -ftop_ksp_type cg -ftop_pc_type jacobi -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
364     filter: grep -v DM_ | grep -v atomic
365 
366 TEST*/
367