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