xref: /petsc/src/dm/impls/swarm/tests/ex2.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 static char help[] = "Tests L2 projection with DMSwarm using delta function particles.\n";
2 
3 #include <petscdmplex.h>
4 #include <petscfe.h>
5 #include <petscdmswarm.h>
6 #include <petscds.h>
7 #include <petscksp.h>
8 #include <petsc/private/petscfeimpl.h>
9 typedef struct {
10   PetscReal L[3];                             /* Dimensions of the mesh bounding box */
11   PetscInt  particlesPerCell;                 /* The number of partices per cell */
12   PetscReal particleRelDx;                    /* Relative particle position perturbation compared to average cell diameter h */
13   PetscReal meshRelDx;                        /* Relative vertex position perturbation compared to average cell diameter h */
14   PetscInt  k;                                /* Mode number for test function */
15   PetscReal momentTol;                        /* Tolerance for checking moment conservation */
16   PetscBool useBlockDiagPrec;                 /* Use the block diagonal of the normal equations as a preconditioner */
17   PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
18 } AppCtx;
19 
20 /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */
21 
22 static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
23 {
24   AppCtx  *ctx = (AppCtx *) a_ctx;
25   PetscInt d;
26 
27   u[0] = 0.0;
28   for (d = 0; d < dim; ++d) u[0] += x[d]/(ctx->L[d]);
29   return 0;
30 }
31 
32 static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
33 {
34   AppCtx  *ctx = (AppCtx *) a_ctx;
35   PetscInt d;
36 
37   u[0] = 1;
38   for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d])*PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4);
39   return 0;
40 }
41 
42 static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
43 {
44   AppCtx *ctx = (AppCtx *) a_ctx;
45 
46   u[0] = PetscSinScalar(2*PETSC_PI*ctx->k*x[0]/(ctx->L[0]));
47   return 0;
48 }
49 
50 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
51 {
52   char           fstring[PETSC_MAX_PATH_LEN] = "linear";
53   PetscBool      flag;
54   PetscErrorCode ierr;
55 
56   PetscFunctionBeginUser;
57   options->particlesPerCell = 1;
58   options->k                = 1;
59   options->particleRelDx    = 1.e-20;
60   options->meshRelDx        = 1.e-20;
61   options->momentTol        = 100.*PETSC_MACHINE_EPSILON;
62   options->useBlockDiagPrec = PETSC_FALSE;
63 
64   ierr = PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");PetscCall(ierr);
65   PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL));
66   PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL));
67   PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL));
68   PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL));
69   PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL));
70   PetscCall(PetscStrcmp(fstring, "linear", &flag));
71   if (flag) {
72     options->func = linear;
73   } else {
74     PetscCall(PetscStrcmp(fstring, "sin", &flag));
75     if (flag) {
76       options->func = sinx;
77     } else {
78       PetscCall(PetscStrcmp(fstring, "x2_x4", &flag));
79       options->func = x2_x4;
80       PetscCheck(flag,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown function %s",fstring);
81     }
82   }
83   PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL));
84   PetscCall(PetscOptionsBool("-block_diag_prec", "Use the block diagonal of the normal equations to precondition the particle projection", "ex2.c", options->useBlockDiagPrec, &options->useBlockDiagPrec, NULL));
85   ierr = PetscOptionsEnd();PetscCall(ierr);
86 
87   PetscFunctionReturn(0);
88 }
89 
90 static PetscErrorCode PerturbVertices(DM dm, AppCtx *user)
91 {
92   PetscRandom    rnd;
93   PetscReal      interval = user->meshRelDx;
94   Vec            coordinates;
95   PetscScalar   *coords;
96   PetscReal     *hh, low[3], high[3];
97   PetscInt       d, cdim, cEnd, N, p, bs;
98 
99   PetscFunctionBeginUser;
100   PetscCall(DMGetBoundingBox(dm, low, high));
101   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
102   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
103   PetscCall(PetscRandomSetInterval(rnd, -interval, interval));
104   PetscCall(PetscRandomSetFromOptions(rnd));
105   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
106   PetscCall(DMGetCoordinateDim(dm, &cdim));
107   PetscCall(PetscCalloc1(cdim,&hh));
108   for (d = 0; d < cdim; ++d) hh[d] = (user->L[d])/PetscPowReal(cEnd, 1./cdim);
109   PetscCall(VecGetLocalSize(coordinates, &N));
110   PetscCall(VecGetBlockSize(coordinates, &bs));
111   PetscCheckFalse(bs != cdim,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %D != %D", bs, cdim);
112   PetscCall(VecGetArray(coordinates, &coords));
113   for (p = 0; p < N; p += cdim) {
114     PetscScalar *coord = &coords[p], value;
115 
116     for (d = 0; d < cdim; ++d) {
117       PetscCall(PetscRandomGetValue(rnd, &value));
118       coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value*hh[d])));
119     }
120   }
121   PetscCall(VecRestoreArray(coordinates, &coords));
122   PetscCall(PetscRandomDestroy(&rnd));
123   PetscCall(PetscFree(hh));
124   PetscFunctionReturn(0);
125 }
126 
127 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
128 {
129   PetscReal      low[3], high[3];
130   PetscInt       cdim, d;
131 
132   PetscFunctionBeginUser;
133   PetscCall(DMCreate(comm, dm));
134   PetscCall(DMSetType(*dm, DMPLEX));
135   PetscCall(DMSetFromOptions(*dm));
136   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
137 
138   PetscCall(DMGetCoordinateDim(*dm, &cdim));
139   PetscCall(DMGetBoundingBox(*dm, low, high));
140   for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d];
141   PetscCall(PerturbVertices(*dm, user));
142   PetscFunctionReturn(0);
143 }
144 
145 static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
146                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
147                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
148                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
149 {
150   g0[0] = 1.0;
151 }
152 
153 static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
154 {
155   PetscFE        fe;
156   PetscDS        ds;
157   DMPolytopeType ct;
158   PetscBool      simplex;
159   PetscInt       dim, cStart;
160 
161   PetscFunctionBeginUser;
162   PetscCall(DMGetDimension(dm, &dim));
163   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
164   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
165   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
166   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, NULL, -1, &fe));
167   PetscCall(PetscObjectSetName((PetscObject) fe, "fe"));
168   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
169   PetscCall(DMCreateDS(dm));
170   PetscCall(PetscFEDestroy(&fe));
171   /* Setup to form mass matrix */
172   PetscCall(DMGetDS(dm, &ds));
173   PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL));
174   PetscFunctionReturn(0);
175 }
176 
177 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
178 {
179   PetscRandom    rnd, rndp;
180   PetscReal      interval = user->particleRelDx;
181   DMPolytopeType ct;
182   PetscBool      simplex;
183   PetscScalar    value, *vals;
184   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
185   PetscInt      *cellid;
186   PetscInt       Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d;
187 
188   PetscFunctionBeginUser;
189   PetscCall(DMGetDimension(dm, &dim));
190   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
191   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
192   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
193   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
194   PetscCall(DMSetType(*sw, DMSWARM));
195   PetscCall(DMSetDimension(*sw, dim));
196 
197   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
198   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
199   PetscCall(PetscRandomSetFromOptions(rnd));
200   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rndp));
201   PetscCall(PetscRandomSetInterval(rndp, -interval, interval));
202   PetscCall(PetscRandomSetFromOptions(rndp));
203 
204   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
205   PetscCall(DMSwarmSetCellDM(*sw, dm));
206   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
207   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
208   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell));
209   PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0));
210   PetscCall(DMSetFromOptions(*sw));
211   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
212   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
213   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **) &vals));
214 
215   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
216   for (c = 0; c < Ncell; ++c) {
217     if (Np == 1) {
218       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
219       cellid[c] = c;
220       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
221     } else {
222       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
223       for (d = 0; d < dim; ++d) xi0[d] = -1.0;
224       for (p = 0; p < Np; ++p) {
225         const PetscInt n   = c*Np + p;
226         PetscReal      sum = 0.0, refcoords[3];
227 
228         cellid[n] = c;
229         for (d = 0; d < dim; ++d) {PetscCall(PetscRandomGetValue(rnd, &value)); refcoords[d] = PetscRealPart(value); sum += refcoords[d];}
230         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
231         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
232       }
233     }
234   }
235   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
236   for (c = 0; c < Ncell; ++c) {
237     for (p = 0; p < Np; ++p) {
238       const PetscInt n = c*Np + p;
239 
240       for (d = 0; d < dim; ++d) {PetscCall(PetscRandomGetValue(rndp, &value)); coords[n*dim+d] += PetscRealPart(value);}
241       user->func(dim, 0.0, &coords[n*dim], 1, &vals[c], user);
242     }
243   }
244   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
245   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
246   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **) &vals));
247   PetscCall(PetscRandomDestroy(&rnd));
248   PetscCall(PetscRandomDestroy(&rndp));
249   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
250   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
251   PetscFunctionReturn(0);
252 }
253 
254 static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user)
255 {
256   DM                 dm;
257   const PetscReal   *coords;
258   const PetscScalar *w;
259   PetscReal          mom[3] = {0.0, 0.0, 0.0};
260   PetscInt           cell, cStart, cEnd, dim;
261 
262   PetscFunctionBeginUser;
263   PetscCall(DMGetDimension(sw, &dim));
264   PetscCall(DMSwarmGetCellDM(sw, &dm));
265   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
266   PetscCall(DMSwarmSortGetAccess(sw));
267   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
268   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &w));
269   for (cell = cStart; cell < cEnd; ++cell) {
270     PetscInt *pidx;
271     PetscInt  Np, p, d;
272 
273     PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
274     for (p = 0; p < Np; ++p) {
275       const PetscInt   idx = pidx[p];
276       const PetscReal *c   = &coords[idx*dim];
277 
278       mom[0] += PetscRealPart(w[idx]);
279       mom[1] += PetscRealPart(w[idx]) * c[0];
280       for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d]*c[d];
281     }
282     PetscCall(PetscFree(pidx));
283   }
284   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
285   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &w));
286   PetscCall(DMSwarmSortRestoreAccess(sw));
287   PetscCallMPI(MPI_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) sw)));
288   PetscFunctionReturn(0);
289 }
290 
291 static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
292                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
293                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
294                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
295 {
296   f0[0] = u[0];
297 }
298 
299 static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux,
300                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
301                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
302                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
303 {
304   f0[0] = x[0]*u[0];
305 }
306 
307 static void f0_r2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
308                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
309                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
310                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
311 {
312   PetscInt d;
313 
314   f0[0] = 0.0;
315   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d])*u[0];
316 }
317 
318 static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
319 {
320   PetscDS        prob;
321   PetscScalar    mom;
322 
323   PetscFunctionBeginUser;
324   PetscCall(DMGetDS(dm, &prob));
325   PetscCall(PetscDSSetObjective(prob, 0, &f0_1));
326   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
327   moments[0] = PetscRealPart(mom);
328   PetscCall(PetscDSSetObjective(prob, 0, &f0_x));
329   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
330   moments[1] = PetscRealPart(mom);
331   PetscCall(PetscDSSetObjective(prob, 0, &f0_r2));
332   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
333   moments[2] = PetscRealPart(mom);
334   PetscFunctionReturn(0);
335 }
336 
337 static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user)
338 {
339   MPI_Comm       comm;
340   KSP            ksp;
341   Mat            M;            /* FEM mass matrix */
342   Mat            M_p;          /* Particle mass matrix */
343   Vec            f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */
344   PetscReal      pmoments[3];  /* \int f, \int x f, \int r^2 f */
345   PetscReal      fmoments[3];  /* \int \hat f, \int x \hat f, \int r^2 \hat f */
346   PetscInt       m;
347 
348   PetscFunctionBeginUser;
349   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
350   PetscCall(KSPCreate(comm, &ksp));
351   PetscCall(KSPSetOptionsPrefix(ksp, "ptof_"));
352   PetscCall(DMGetGlobalVector(dm, &fhat));
353   PetscCall(DMGetGlobalVector(dm, &rhs));
354 
355   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
356   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
357 
358   /* make particle weight vector */
359   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
360 
361   /* create matrix RHS vector */
362   PetscCall(MatMultTranspose(M_p, f, rhs));
363   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
364   PetscCall(PetscObjectSetName((PetscObject) rhs,"rhs"));
365   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
366 
367   PetscCall(DMCreateMatrix(dm, &M));
368   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
369   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
370   PetscCall(KSPSetOperators(ksp, M, M));
371   PetscCall(KSPSetFromOptions(ksp));
372   PetscCall(KSPSolve(ksp, rhs, fhat));
373   PetscCall(PetscObjectSetName((PetscObject) fhat,"fhat"));
374   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
375 
376   /* Check moments of field */
377   PetscCall(computeParticleMoments(sw, pmoments, user));
378   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
379   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
380   for (m = 0; m < 3; ++m) {
381     PetscCheckFalse(PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol,comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol);
382   }
383 
384   PetscCall(KSPDestroy(&ksp));
385   PetscCall(MatDestroy(&M));
386   PetscCall(MatDestroy(&M_p));
387   PetscCall(DMRestoreGlobalVector(dm, &fhat));
388   PetscCall(DMRestoreGlobalVector(dm, &rhs));
389 
390   PetscFunctionReturn(0);
391 }
392 
393 static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user)
394 {
395 
396   MPI_Comm       comm;
397   KSP            ksp;
398   Mat            M;            /* FEM mass matrix */
399   Mat            M_p, PM_p;    /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */
400   Vec            f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */
401   PetscReal      pmoments[3];  /* \int f, \int x f, \int r^2 f */
402   PetscReal      fmoments[3];  /* \int \hat f, \int x \hat f, \int r^2 \hat f */
403   PetscInt       m;
404 
405   PetscFunctionBeginUser;
406   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
407 
408   PetscCall(KSPCreate(comm, &ksp));
409   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
410   PetscCall(KSPSetFromOptions(ksp));
411 
412   PetscCall(DMGetGlobalVector(dm, &fhat));
413   PetscCall(DMGetGlobalVector(dm, &rhs));
414 
415   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
416   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
417 
418   /* make particle weight vector */
419   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
420 
421   /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */
422   PetscCall(PetscObjectSetName((PetscObject) rhs,"rhs"));
423   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
424   PetscCall(DMCreateMatrix(dm, &M));
425   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
426   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
427   PetscCall(MatMultTranspose(M, fhat, rhs));
428   if (user->useBlockDiagPrec) PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
429   else                        {PetscCall(PetscObjectReference((PetscObject) M_p)); PM_p = M_p;}
430 
431   PetscCall(KSPSetOperators(ksp, M_p, PM_p));
432   PetscCall(KSPSolveTranspose(ksp, rhs, f));
433   PetscCall(PetscObjectSetName((PetscObject) fhat,"fhat"));
434   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
435 
436   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
437 
438   /* Check moments */
439   PetscCall(computeParticleMoments(sw, pmoments, user));
440   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
441   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
442   for (m = 0; m < 3; ++m) {
443     PetscCheckFalse(PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol,comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol);
444   }
445   PetscCall(KSPDestroy(&ksp));
446   PetscCall(MatDestroy(&M));
447   PetscCall(MatDestroy(&M_p));
448   PetscCall(MatDestroy(&PM_p));
449   PetscCall(DMRestoreGlobalVector(dm, &fhat));
450   PetscCall(DMRestoreGlobalVector(dm, &rhs));
451   PetscFunctionReturn(0);
452 }
453 
454 /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */
455 static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC)
456 {
457   DM_Plex         *mesh  = (DM_Plex *) dm->data;
458   PetscInt         debug = mesh->printFEM;
459   DM               dmC;
460   PetscSection     section;
461   PetscQuadrature  quad = NULL;
462   PetscScalar     *interpolant, *gradsum;
463   PetscFEGeom      fegeom;
464   PetscReal       *coords;
465   const PetscReal *quadPoints, *quadWeights;
466   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;
467 
468   PetscFunctionBegin;
469   PetscCall(VecGetDM(locC, &dmC));
470   PetscCall(VecSet(locC, 0.0));
471   PetscCall(DMGetDimension(dm, &dim));
472   PetscCall(DMGetCoordinateDim(dm, &coordDim));
473   fegeom.dimEmbed = coordDim;
474   PetscCall(DMGetLocalSection(dm, &section));
475   PetscCall(PetscSectionGetNumFields(section, &numFields));
476   for (field = 0; field < numFields; ++field) {
477     PetscObject  obj;
478     PetscClassId id;
479     PetscInt     Nc;
480 
481     PetscCall(DMGetField(dm, field, NULL, &obj));
482     PetscCall(PetscObjectGetClassId(obj, &id));
483     if (id == PETSCFE_CLASSID) {
484       PetscFE fe = (PetscFE) obj;
485 
486       PetscCall(PetscFEGetQuadrature(fe, &quad));
487       PetscCall(PetscFEGetNumComponents(fe, &Nc));
488     } else if (id == PETSCFV_CLASSID) {
489       PetscFV fv = (PetscFV) obj;
490 
491       PetscCall(PetscFVGetQuadrature(fv, &quad));
492       PetscCall(PetscFVGetNumComponents(fv, &Nc));
493     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
494     numComponents += Nc;
495   }
496   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
497   PetscCheckFalse((qNc != 1) && (qNc != numComponents),PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
498   PetscCall(PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ));
499   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
500   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
501   for (v = vStart; v < vEnd; ++v) {
502     PetscInt *star = NULL;
503     PetscInt starSize, st, d, fc;
504 
505     PetscCall(PetscArrayzero(gradsum, coordDim*numComponents));
506     PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
507     for (st = 0; st < starSize*2; st += 2) {
508       const PetscInt cell = star[st];
509       PetscScalar    *grad = &gradsum[coordDim*numComponents];
510       PetscScalar    *x    = NULL;
511 
512       if ((cell < cStart) || (cell >= cEnd)) continue;
513       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
514       PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
515       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
516         PetscObject  obj;
517         PetscClassId id;
518         PetscInt     Nb, Nc, q, qc = 0;
519 
520         PetscCall(PetscArrayzero(grad, coordDim*numComponents));
521         PetscCall(DMGetField(dm, field, NULL, &obj));
522         PetscCall(PetscObjectGetClassId(obj, &id));
523         if (id == PETSCFE_CLASSID)      {PetscCall(PetscFEGetNumComponents((PetscFE) obj, &Nc));PetscCall(PetscFEGetDimension((PetscFE) obj, &Nb));}
524         else if (id == PETSCFV_CLASSID) {PetscCall(PetscFVGetNumComponents((PetscFV) obj, &Nc));Nb = 1;}
525         else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
526         for (q = 0; q < Nq; ++q) {
527           PetscCheck(fegeom.detJ[q] > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)fegeom.detJ[q], cell, q);
528           if (id == PETSCFE_CLASSID)      PetscCall(PetscFEInterpolateGradient_Static((PetscFE) obj, 1, &x[fieldOffset], &fegeom, q, interpolant));
529           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
530           for (fc = 0; fc < Nc; ++fc) {
531             const PetscReal wt = quadWeights[q*qNc+qc+fc];
532 
533             for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*fegeom.detJ[q];
534           }
535         }
536         fieldOffset += Nb;
537       }
538       PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
539       for (fc = 0; fc < numComponents; ++fc) {
540         for (d = 0; d < coordDim; ++d) {
541           gradsum[fc*coordDim+d] += grad[fc*coordDim+d];
542         }
543       }
544       if (debug) {
545         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell));
546         for (fc = 0; fc < numComponents; ++fc) {
547           for (d = 0; d < coordDim; ++d) {
548             if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
549             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d])));
550           }
551         }
552         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
553       }
554     }
555     PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
556     PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
557   }
558   PetscCall(PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ));
559   PetscFunctionReturn(0);
560 }
561 
562 static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user)
563 {
564 
565   MPI_Comm       comm;
566   KSP            ksp;
567   Mat            M;                   /* FEM mass matrix */
568   Mat            M_p;                 /* Particle mass matrix */
569   Vec            f, rhs, fhat, grad;  /* Particle field f, \int phi_i f, FEM field */
570   PetscReal      pmoments[3];         /* \int f, \int x f, \int r^2 f */
571   PetscReal      fmoments[3];         /* \int \hat f, \int x \hat f, \int r^2 \hat f */
572   PetscInt       m;
573 
574   PetscFunctionBeginUser;
575   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
576   PetscCall(KSPCreate(comm, &ksp));
577   PetscCall(KSPSetOptionsPrefix(ksp, "ptof_"));
578   PetscCall(DMGetGlobalVector(dm, &fhat));
579   PetscCall(DMGetGlobalVector(dm, &rhs));
580   PetscCall(DMGetGlobalVector(dm, &grad));
581 
582   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
583   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
584 
585   /* make particle weight vector */
586   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
587 
588   /* create matrix RHS vector */
589   PetscCall(MatMultTranspose(M_p, f, rhs));
590   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
591   PetscCall(PetscObjectSetName((PetscObject) rhs,"rhs"));
592   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
593 
594   PetscCall(DMCreateMatrix(dm, &M));
595   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
596 
597   PetscCall(InterpolateGradient(dm, fhat, grad));
598 
599   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
600   PetscCall(KSPSetOperators(ksp, M, M));
601   PetscCall(KSPSetFromOptions(ksp));
602   PetscCall(KSPSolve(ksp, rhs, grad));
603   PetscCall(PetscObjectSetName((PetscObject) fhat,"fhat"));
604   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
605 
606   /* Check moments of field */
607   PetscCall(computeParticleMoments(sw, pmoments, user));
608   PetscCall(computeFEMMoments(dm, grad, fmoments, user));
609   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
610   for (m = 0; m < 3; ++m) {
611     PetscCheckFalse(PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol,comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol);
612   }
613 
614   PetscCall(KSPDestroy(&ksp));
615   PetscCall(MatDestroy(&M));
616   PetscCall(MatDestroy(&M_p));
617   PetscCall(DMRestoreGlobalVector(dm, &fhat));
618   PetscCall(DMRestoreGlobalVector(dm, &rhs));
619   PetscCall(DMRestoreGlobalVector(dm, &grad));
620 
621   PetscFunctionReturn(0);
622 }
623 
624 int main (int argc, char *argv[])
625 {
626   MPI_Comm       comm;
627   DM             dm, sw;
628   AppCtx         user;
629 
630   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
631   comm = PETSC_COMM_WORLD;
632   PetscCall(ProcessOptions(comm, &user));
633   PetscCall(CreateMesh(comm, &dm, &user));
634   PetscCall(CreateFEM(dm, &user));
635   PetscCall(CreateParticles(dm, &sw, &user));
636   PetscCall(TestL2ProjectionParticlesToField(dm, sw, &user));
637   PetscCall(TestL2ProjectionFieldToParticles(dm, sw, &user));
638   PetscCall(TestFieldGradientProjection(dm, sw, &user));
639   PetscCall(DMDestroy(&dm));
640   PetscCall(DMDestroy(&sw));
641   PetscCall(PetscFinalize());
642   return 0;
643 }
644 
645 /*TEST
646 
647   # Swarm does not handle complex or quad
648   build:
649     requires: !complex double
650 
651   test:
652     suffix: proj_tri_0
653     requires: triangle
654     args: -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
655     filter: grep -v marker | grep -v atomic | grep -v usage
656 
657   test:
658     suffix: proj_tri_2_faces
659     requires: triangle
660     args: -dm_plex_box_faces 2,2  -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
661     filter: grep -v marker | grep -v atomic | grep -v usage
662 
663   test:
664     suffix: proj_quad_0
665     requires: triangle
666     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
667     filter: grep -v marker | grep -v atomic | grep -v usage
668 
669   test:
670     suffix: proj_quad_2_faces
671     requires: triangle
672     args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
673     filter: grep -v marker | grep -v atomic | grep -v usage
674 
675   test:
676     suffix: proj_tri_5P
677     requires: triangle
678     args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
679     filter: grep -v marker | grep -v atomic | grep -v usage
680 
681   test:
682     suffix: proj_quad_5P
683     requires: triangle
684     args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
685     filter: grep -v marker | grep -v atomic | grep -v usage
686 
687   test:
688     suffix: proj_tri_mdx
689     requires: triangle
690     args: -dm_plex_box_faces 1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
691     filter: grep -v marker | grep -v atomic | grep -v usage
692 
693   test:
694     suffix: proj_tri_mdx_5P
695     requires: triangle
696     args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
697     filter: grep -v marker | grep -v atomic | grep -v usage
698 
699   test:
700     suffix: proj_tri_3d
701     requires: ctetgen
702     args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
703     filter: grep -v marker | grep -v atomic | grep -v usage
704 
705   test:
706     suffix: proj_tri_3d_2_faces
707     requires: ctetgen
708     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
709     filter: grep -v marker | grep -v atomic | grep -v usage
710 
711   test:
712     suffix: proj_tri_3d_5P
713     requires: ctetgen
714     args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
715     filter: grep -v marker | grep -v atomic | grep -v usage
716 
717   test:
718     suffix: proj_tri_3d_mdx
719     requires: ctetgen
720     args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
721     filter: grep -v marker | grep -v atomic | grep -v usage
722 
723   test:
724     suffix: proj_tri_3d_mdx_5P
725     requires: ctetgen
726     args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
727     filter: grep -v marker | grep -v atomic | grep -v usage
728 
729   test:
730     suffix: proj_tri_3d_mdx_2_faces
731     requires: ctetgen
732     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
733     filter: grep -v marker | grep -v atomic | grep -v usage
734 
735   test:
736     suffix: proj_tri_3d_mdx_5P_2_faces
737     requires: ctetgen
738     args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu  -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
739     filter: grep -v marker | grep -v atomic | grep -v usage
740 
741   test:
742     suffix: proj_quad_lsqr_scale
743     requires: !complex
744     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
745       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
746       -particlesPerCell 200 \
747       -ptof_pc_type lu  \
748       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none
749 
750   test:
751     suffix: proj_quad_lsqr_prec_scale
752     requires: !complex
753     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
754       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
755       -particlesPerCell 200 \
756       -ptof_pc_type lu  \
757       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec
758 
759 TEST*/
760