xref: /petsc/src/dm/impls/swarm/tests/ex2.c (revision bcee047adeeb73090d7e36cc71e39fc287cdbb97)
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 useBlockDiagPrec; /* Use the block diagonal of the normal equations as a preconditioner */
17   PetscBool shape;
18   PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
19 } AppCtx;
20 
21 /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */
22 
23 static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
24 {
25   AppCtx  *ctx = (AppCtx *)a_ctx;
26   PetscInt d;
27 
28   u[0] = 0.0;
29   for (d = 0; d < dim; ++d) u[0] += x[d] / (ctx->L[d]);
30   return PETSC_SUCCESS;
31 }
32 
33 static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
34 {
35   AppCtx  *ctx = (AppCtx *)a_ctx;
36   PetscInt d;
37 
38   u[0] = 1;
39   for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d]) * PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4);
40   return PETSC_SUCCESS;
41 }
42 
43 static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
44 {
45   AppCtx *ctx = (AppCtx *)a_ctx;
46 
47   u[0] = PetscSinScalar(2 * PETSC_PI * ctx->k * x[0] / (ctx->L[0]));
48   return PETSC_SUCCESS;
49 }
50 
51 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
52 {
53   char      fstring[PETSC_MAX_PATH_LEN] = "linear";
54   PetscBool flag;
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   options->shape            = PETSC_FALSE;
64 
65   PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");
66   PetscCall(PetscOptionsBool("-shape", "Test linear function projection", "ex2.c", options->shape, &options->shape, NULL));
67   PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL));
68   PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL));
69   PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL));
70   PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL));
71   PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL));
72   PetscCall(PetscStrcmp(fstring, "linear", &flag));
73   if (flag) {
74     options->func = linear;
75   } else {
76     PetscCall(PetscStrcmp(fstring, "sin", &flag));
77     if (flag) {
78       options->func = sinx;
79     } else {
80       PetscCall(PetscStrcmp(fstring, "x2_x4", &flag));
81       options->func = x2_x4;
82       PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown function %s", fstring);
83     }
84   }
85   PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL));
86   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));
87   PetscOptionsEnd();
88 
89   PetscFunctionReturn(PETSC_SUCCESS);
90 }
91 
92 static PetscErrorCode PerturbVertices(DM dm, AppCtx *user)
93 {
94   PetscRandom  rnd;
95   PetscReal    interval = user->meshRelDx;
96   Vec          coordinates;
97   PetscScalar *coords;
98   PetscReal   *hh, low[3], high[3];
99   PetscInt     d, cdim, cEnd, N, p, bs;
100 
101   PetscFunctionBeginUser;
102   PetscCall(DMGetBoundingBox(dm, low, high));
103   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
104   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
105   PetscCall(PetscRandomSetInterval(rnd, -interval, interval));
106   PetscCall(PetscRandomSetFromOptions(rnd));
107   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
108   PetscCall(DMGetCoordinateDim(dm, &cdim));
109   PetscCall(PetscCalloc1(cdim, &hh));
110   for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim);
111   PetscCall(VecGetLocalSize(coordinates, &N));
112   PetscCall(VecGetBlockSize(coordinates, &bs));
113   PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim);
114   PetscCall(VecGetArray(coordinates, &coords));
115   for (p = 0; p < N; p += cdim) {
116     PetscScalar *coord = &coords[p], value;
117 
118     for (d = 0; d < cdim; ++d) {
119       PetscCall(PetscRandomGetValue(rnd, &value));
120       coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d])));
121     }
122   }
123   PetscCall(VecRestoreArray(coordinates, &coords));
124   PetscCall(PetscRandomDestroy(&rnd));
125   PetscCall(PetscFree(hh));
126   PetscFunctionReturn(PETSC_SUCCESS);
127 }
128 
129 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
130 {
131   PetscReal low[3], high[3];
132   PetscInt  cdim, d;
133 
134   PetscFunctionBeginUser;
135   PetscCall(DMCreate(comm, dm));
136   PetscCall(DMSetType(*dm, DMPLEX));
137   PetscCall(DMSetFromOptions(*dm));
138   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
139 
140   PetscCall(DMGetCoordinateDim(*dm, &cdim));
141   PetscCall(DMGetBoundingBox(*dm, low, high));
142   for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d];
143   PetscCall(PerturbVertices(*dm, user));
144   PetscFunctionReturn(PETSC_SUCCESS);
145 }
146 
147 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[])
148 {
149   g0[0] = 1.0;
150 }
151 
152 static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
153 {
154   PetscFE        fe;
155   PetscDS        ds;
156   DMPolytopeType ct;
157   PetscBool      simplex;
158   PetscInt       dim, cStart;
159 
160   PetscFunctionBeginUser;
161   PetscCall(DMGetDimension(dm, &dim));
162   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
163   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
164   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
165   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe));
166   PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
167   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
168   PetscCall(DMCreateDS(dm));
169   PetscCall(PetscFEDestroy(&fe));
170   /* Setup to form mass matrix */
171   PetscCall(DMGetDS(dm, &ds));
172   PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL));
173   PetscFunctionReturn(PETSC_SUCCESS);
174 }
175 
176 static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
177 {
178   PetscRandom    rnd, rndp;
179   PetscReal      interval = user->particleRelDx;
180   DMPolytopeType ct;
181   PetscBool      simplex;
182   PetscScalar    value, *vals;
183   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
184   PetscInt      *cellid;
185   PetscInt       Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d;
186 
187   PetscFunctionBeginUser;
188   PetscCall(DMGetDimension(dm, &dim));
189   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
190   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
191   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
192   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
193   PetscCall(DMSetType(*sw, DMSWARM));
194   PetscCall(DMSetDimension(*sw, dim));
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, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
211   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&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       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         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, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
252   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&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(DMViewFromOptions(*sw, NULL, "-sw_view"));
258   PetscFunctionReturn(PETSC_SUCCESS);
259 }
260 
261 static PetscErrorCode CreateParticles_Shape(DM dm, DM *sw, AppCtx *user)
262 {
263   PetscDS          prob;
264   PetscFE          fe;
265   PetscQuadrature  quad;
266   PetscScalar     *vals;
267   PetscReal       *v0, *J, *invJ, detJ, *coords, *xi0;
268   PetscInt        *cellid;
269   const PetscReal *qpoints, *qweights;
270   PetscInt         cStart, cEnd, c, Nq, q, dim;
271 
272   PetscFunctionBeginUser;
273   PetscCall(DMGetDimension(dm, &dim));
274   PetscCall(DMGetDS(dm, &prob));
275   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
276   PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
277   PetscCall(PetscFEGetQuadrature(fe, &quad));
278   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, &qweights));
279   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
280   PetscCall(DMSetType(*sw, DMSWARM));
281   PetscCall(DMSetDimension(*sw, dim));
282   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
283   PetscCall(DMSwarmSetCellDM(*sw, dm));
284   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_REAL));
285   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
286   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Nq, 0));
287   PetscCall(DMSetFromOptions(*sw));
288   PetscCall(PetscMalloc4(dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
289   for (c = 0; c < dim; c++) xi0[c] = -1.;
290   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
291   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
292   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals));
293   for (c = cStart; c < cEnd; ++c) {
294     for (q = 0; q < Nq; ++q) {
295       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
296       cellid[c * Nq + q] = c;
297       CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], &coords[(c * Nq + q) * dim]);
298       PetscCall(linear(dim, 0.0, &coords[(c * Nq + q) * dim], 1, &vals[c * Nq + q], (void *)user));
299       vals[c * Nq + q] *= qweights[q] * detJ;
300     }
301   }
302   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
303   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
304   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals));
305   PetscCall(PetscFree4(xi0, v0, J, invJ));
306   PetscCall(DMSwarmMigrate(*sw, PETSC_FALSE));
307   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
308   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
309   PetscFunctionReturn(PETSC_SUCCESS);
310 }
311 
312 static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user)
313 {
314   DM                 dm;
315   const PetscReal   *coords;
316   const PetscScalar *w;
317   PetscReal          mom[3] = {0.0, 0.0, 0.0};
318   PetscInt           cell, cStart, cEnd, dim;
319 
320   PetscFunctionBeginUser;
321   PetscCall(DMGetDimension(sw, &dim));
322   PetscCall(DMSwarmGetCellDM(sw, &dm));
323   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
324   PetscCall(DMSwarmSortGetAccess(sw));
325   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
326   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w));
327   for (cell = cStart; cell < cEnd; ++cell) {
328     PetscInt *pidx;
329     PetscInt  Np, p, d;
330 
331     PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
332     for (p = 0; p < Np; ++p) {
333       const PetscInt   idx = pidx[p];
334       const PetscReal *c   = &coords[idx * dim];
335 
336       mom[0] += PetscRealPart(w[idx]);
337       mom[1] += PetscRealPart(w[idx]) * c[0];
338       for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d];
339     }
340     PetscCall(PetscFree(pidx));
341   }
342   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
343   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
344   PetscCall(DMSwarmSortRestoreAccess(sw));
345   PetscCall(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
346   PetscFunctionReturn(PETSC_SUCCESS);
347 }
348 
349 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[])
350 {
351   f0[0] = u[0];
352 }
353 
354 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[])
355 {
356   f0[0] = x[0] * u[0];
357 }
358 
359 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[])
360 {
361   PetscInt d;
362 
363   f0[0] = 0.0;
364   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0];
365 }
366 
367 static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
368 {
369   PetscDS     prob;
370   PetscScalar mom;
371 
372   PetscFunctionBeginUser;
373   PetscCall(DMGetDS(dm, &prob));
374   PetscCall(PetscDSSetObjective(prob, 0, &f0_1));
375   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
376   moments[0] = PetscRealPart(mom);
377   PetscCall(PetscDSSetObjective(prob, 0, &f0_x));
378   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
379   moments[1] = PetscRealPart(mom);
380   PetscCall(PetscDSSetObjective(prob, 0, &f0_r2));
381   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
382   moments[2] = PetscRealPart(mom);
383   PetscFunctionReturn(PETSC_SUCCESS);
384 }
385 
386 static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user)
387 {
388   MPI_Comm  comm;
389   KSP       ksp;
390   Mat       M;            /* FEM mass matrix */
391   Mat       M_p;          /* Particle mass matrix */
392   Vec       f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */
393   PetscReal pmoments[3];  /* \int f, \int x f, \int r^2 f */
394   PetscReal fmoments[3];  /* \int \hat f, \int x \hat f, \int r^2 \hat f */
395   PetscInt  m;
396 
397   PetscFunctionBeginUser;
398   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
399   PetscCall(KSPCreate(comm, &ksp));
400   PetscCall(KSPSetOptionsPrefix(ksp, "ptof_"));
401   PetscCall(DMGetGlobalVector(dm, &fhat));
402   PetscCall(DMGetGlobalVector(dm, &rhs));
403 
404   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
405   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
406 
407   /* make particle weight vector */
408   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
409 
410   /* create matrix RHS vector */
411   PetscCall(MatMultTranspose(M_p, f, rhs));
412   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
413   PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs"));
414   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
415 
416   PetscCall(DMCreateMatrix(dm, &M));
417   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
418   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
419   PetscCall(KSPSetOperators(ksp, M, M));
420   PetscCall(KSPSetFromOptions(ksp));
421   PetscCall(KSPSolve(ksp, rhs, fhat));
422   PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat"));
423   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
424 
425   /* Check moments of field */
426   PetscCall(computeParticleMoments(sw, pmoments, user));
427   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
428   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
429   for (m = 0; m < 3; ++m) {
430     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);
431   }
432 
433   PetscCall(KSPDestroy(&ksp));
434   PetscCall(MatDestroy(&M));
435   PetscCall(MatDestroy(&M_p));
436   PetscCall(DMRestoreGlobalVector(dm, &fhat));
437   PetscCall(DMRestoreGlobalVector(dm, &rhs));
438 
439   PetscFunctionReturn(PETSC_SUCCESS);
440 }
441 
442 static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user)
443 {
444   MPI_Comm  comm;
445   KSP       ksp;
446   Mat       M;            /* FEM mass matrix */
447   Mat       M_p, PM_p;    /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */
448   Vec       f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */
449   PetscReal pmoments[3];  /* \int f, \int x f, \int r^2 f */
450   PetscReal fmoments[3];  /* \int \hat f, \int x \hat f, \int r^2 \hat f */
451   PetscInt  m;
452 
453   PetscFunctionBeginUser;
454   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
455 
456   PetscCall(KSPCreate(comm, &ksp));
457   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
458   PetscCall(KSPSetFromOptions(ksp));
459   PetscCall(DMGetGlobalVector(dm, &fhat));
460   PetscCall(DMGetGlobalVector(dm, &rhs));
461   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
462   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
463   /* make particle weight vector */
464   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
465   /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */
466   PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs"));
467   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
468   PetscCall(DMCreateMatrix(dm, &M));
469   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
470   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
471   PetscCall(MatMultTranspose(M, fhat, rhs));
472   if (user->useBlockDiagPrec) PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
473   else {
474     PetscCall(PetscObjectReference((PetscObject)M_p));
475     PM_p = M_p;
476   }
477   PetscCall(KSPSetOperators(ksp, M_p, PM_p));
478   PetscCall(KSPSolveTranspose(ksp, rhs, f));
479   PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat"));
480   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
481   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
482   /* Check moments */
483   PetscCall(computeParticleMoments(sw, pmoments, user));
484   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
485   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
486   for (m = 0; m < 3; ++m) {
487     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);
488   }
489   PetscCall(KSPDestroy(&ksp));
490   PetscCall(MatDestroy(&M));
491   PetscCall(MatDestroy(&M_p));
492   PetscCall(MatDestroy(&PM_p));
493   PetscCall(DMRestoreGlobalVector(dm, &fhat));
494   PetscCall(DMRestoreGlobalVector(dm, &rhs));
495   PetscFunctionReturn(PETSC_SUCCESS);
496 }
497 
498 /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */
499 static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC)
500 {
501   DM_Plex         *mesh  = (DM_Plex *)dm->data;
502   PetscInt         debug = mesh->printFEM;
503   DM               dmC;
504   PetscSection     section;
505   PetscQuadrature  quad = NULL;
506   PetscScalar     *interpolant, *gradsum;
507   PetscFEGeom      fegeom;
508   PetscReal       *coords;
509   const PetscReal *quadPoints, *quadWeights;
510   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;
511 
512   PetscFunctionBegin;
513   PetscCall(VecGetDM(locC, &dmC));
514   PetscCall(VecSet(locC, 0.0));
515   PetscCall(DMGetDimension(dm, &dim));
516   PetscCall(DMGetCoordinateDim(dm, &coordDim));
517   fegeom.dimEmbed = coordDim;
518   PetscCall(DMGetLocalSection(dm, &section));
519   PetscCall(PetscSectionGetNumFields(section, &numFields));
520   for (field = 0; field < numFields; ++field) {
521     PetscObject  obj;
522     PetscClassId id;
523     PetscInt     Nc;
524 
525     PetscCall(DMGetField(dm, field, NULL, &obj));
526     PetscCall(PetscObjectGetClassId(obj, &id));
527     if (id == PETSCFE_CLASSID) {
528       PetscFE fe = (PetscFE)obj;
529 
530       PetscCall(PetscFEGetQuadrature(fe, &quad));
531       PetscCall(PetscFEGetNumComponents(fe, &Nc));
532     } else if (id == PETSCFV_CLASSID) {
533       PetscFV fv = (PetscFV)obj;
534 
535       PetscCall(PetscFVGetQuadrature(fv, &quad));
536       PetscCall(PetscFVGetNumComponents(fv, &Nc));
537     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
538     numComponents += Nc;
539   }
540   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
541   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
542   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));
543   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
544   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
545   for (v = vStart; v < vEnd; ++v) {
546     PetscInt *star = NULL;
547     PetscInt  starSize, st, d, fc;
548 
549     PetscCall(PetscArrayzero(gradsum, coordDim * numComponents));
550     PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
551     for (st = 0; st < starSize * 2; st += 2) {
552       const PetscInt cell = star[st];
553       PetscScalar   *grad = &gradsum[coordDim * numComponents];
554       PetscScalar   *x    = NULL;
555 
556       if ((cell < cStart) || (cell >= cEnd)) continue;
557       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
558       PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
559       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
560         PetscObject  obj;
561         PetscClassId id;
562         PetscInt     Nb, Nc, q, qc = 0;
563 
564         PetscCall(PetscArrayzero(grad, coordDim * numComponents));
565         PetscCall(DMGetField(dm, field, NULL, &obj));
566         PetscCall(PetscObjectGetClassId(obj, &id));
567         if (id == PETSCFE_CLASSID) {
568           PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
569           PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
570         } else if (id == PETSCFV_CLASSID) {
571           PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
572           Nb = 1;
573         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
574         for (q = 0; q < Nq; ++q) {
575           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);
576           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant));
577           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
578           for (fc = 0; fc < Nc; ++fc) {
579             const PetscReal wt = quadWeights[q * qNc + qc + fc];
580 
581             for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q];
582           }
583         }
584         fieldOffset += Nb;
585       }
586       PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
587       for (fc = 0; fc < numComponents; ++fc) {
588         for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d];
589       }
590       if (debug) {
591         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell));
592         for (fc = 0; fc < numComponents; ++fc) {
593           for (d = 0; d < coordDim; ++d) {
594             if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
595             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d])));
596           }
597         }
598         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
599       }
600     }
601     PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
602     PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
603   }
604   PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
605   PetscFunctionReturn(PETSC_SUCCESS);
606 }
607 
608 static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user)
609 {
610   MPI_Comm  comm;
611   KSP       ksp;
612   Mat       M;                  /* FEM mass matrix */
613   Mat       M_p;                /* Particle mass matrix */
614   Vec       f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */
615   PetscReal pmoments[3];        /* \int f, \int x f, \int r^2 f */
616   PetscReal fmoments[3];        /* \int \hat f, \int x \hat f, \int r^2 \hat f */
617   PetscInt  m;
618 
619   PetscFunctionBeginUser;
620   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
621   PetscCall(KSPCreate(comm, &ksp));
622   PetscCall(KSPSetOptionsPrefix(ksp, "ptof_"));
623   PetscCall(DMGetGlobalVector(dm, &fhat));
624   PetscCall(DMGetGlobalVector(dm, &rhs));
625   PetscCall(DMGetGlobalVector(dm, &grad));
626 
627   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
628   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
629 
630   /* make particle weight vector */
631   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
632 
633   /* create matrix RHS vector */
634   PetscCall(MatMultTranspose(M_p, f, rhs));
635   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
636   PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs"));
637   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
638 
639   PetscCall(DMCreateMatrix(dm, &M));
640   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
641 
642   PetscCall(InterpolateGradient(dm, fhat, grad));
643 
644   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
645   PetscCall(KSPSetOperators(ksp, M, M));
646   PetscCall(KSPSetFromOptions(ksp));
647   PetscCall(KSPSolve(ksp, rhs, grad));
648   PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat"));
649   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
650 
651   /* Check moments of field */
652   PetscCall(computeParticleMoments(sw, pmoments, user));
653   PetscCall(computeFEMMoments(dm, grad, fmoments, user));
654   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
655   for (m = 0; m < 3; ++m) {
656     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);
657   }
658 
659   PetscCall(KSPDestroy(&ksp));
660   PetscCall(MatDestroy(&M));
661   PetscCall(MatDestroy(&M_p));
662   PetscCall(DMRestoreGlobalVector(dm, &fhat));
663   PetscCall(DMRestoreGlobalVector(dm, &rhs));
664   PetscCall(DMRestoreGlobalVector(dm, &grad));
665 
666   PetscFunctionReturn(PETSC_SUCCESS);
667 }
668 
669 static PetscErrorCode TestL2Projection_Shape(DM dm, DM sw, AppCtx *user)
670 {
671   PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
672   KSP       ksp;
673   Mat       mass;
674   Vec       u, rhs, uproj;
675   void     *ctxs[1];
676   PetscReal error;
677 
678   PetscFunctionBeginUser;
679   ctxs[0]  = (void *)user;
680   funcs[0] = linear;
681 
682   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &u));
683   PetscCall(VecViewFromOptions(u, NULL, "-f_view"));
684   PetscCall(DMGetGlobalVector(dm, &rhs));
685   PetscCall(DMCreateMassMatrix(sw, dm, &mass));
686   PetscCall(MatMultTranspose(mass, u, rhs));
687   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
688   PetscCall(MatDestroy(&mass));
689   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &u));
690   PetscCall(DMGetGlobalVector(dm, &uproj));
691   PetscCall(DMCreateMatrix(dm, &mass));
692   PetscCall(DMPlexSNESComputeJacobianFEM(dm, uproj, mass, mass, user));
693   PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view"));
694   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
695   PetscCall(KSPSetOperators(ksp, mass, mass));
696   PetscCall(KSPSetFromOptions(ksp));
697   PetscCall(KSPSolve(ksp, rhs, uproj));
698   PetscCall(KSPDestroy(&ksp));
699   PetscCall(PetscObjectSetName((PetscObject)uproj, "Full Projection"));
700   PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view"));
701   PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, uproj, &error));
702   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error));
703   PetscCall(MatDestroy(&mass));
704   PetscCall(DMRestoreGlobalVector(dm, &rhs));
705   PetscCall(DMRestoreGlobalVector(dm, &uproj));
706   PetscFunctionReturn(PETSC_SUCCESS);
707 }
708 
709 int main(int argc, char *argv[])
710 {
711   MPI_Comm comm;
712   DM       dm, sw;
713   AppCtx   user;
714 
715   PetscFunctionBeginUser;
716   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
717   comm = PETSC_COMM_WORLD;
718   PetscCall(ProcessOptions(comm, &user));
719   PetscCall(CreateMesh(comm, &dm, &user));
720   PetscCall(CreateFEM(dm, &user));
721   if (!user.shape) {
722     PetscCall(CreateParticles(dm, &sw, &user));
723     PetscCall(TestL2ProjectionParticlesToField(dm, sw, &user));
724     PetscCall(TestL2ProjectionFieldToParticles(dm, sw, &user));
725     PetscCall(TestFieldGradientProjection(dm, sw, &user));
726   } else {
727     PetscCall(CreateParticles_Shape(dm, &sw, &user));
728     PetscCall(TestL2Projection_Shape(dm, sw, &user));
729   }
730   PetscCall(DMDestroy(&dm));
731   PetscCall(DMDestroy(&sw));
732   PetscCall(PetscFinalize());
733   return 0;
734 }
735 
736 /*TEST
737 
738   # Swarm does not handle complex or quad
739   build:
740     requires: !complex double
741 
742   test:
743     suffix: proj_tri_0
744     requires: triangle
745     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
746     filter: grep -v marker | grep -v atomic | grep -v usage
747 
748   test:
749     suffix: proj_tri_2_faces
750     requires: triangle
751     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
752     filter: grep -v marker | grep -v atomic | grep -v usage
753 
754   test:
755     suffix: proj_quad_0
756     requires: triangle
757     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
758     filter: grep -v marker | grep -v atomic | grep -v usage
759 
760   test:
761     suffix: proj_quad_2_faces
762     requires: triangle
763     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
764     filter: grep -v marker | grep -v atomic | grep -v usage
765 
766   test:
767     suffix: proj_tri_5P
768     requires: triangle
769     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
770     filter: grep -v marker | grep -v atomic | grep -v usage
771 
772   test:
773     suffix: proj_quad_5P
774     requires: triangle
775     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
776     filter: grep -v marker | grep -v atomic | grep -v usage
777 
778   test:
779     suffix: proj_tri_mdx
780     requires: triangle
781     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
782     filter: grep -v marker | grep -v atomic | grep -v usage
783 
784   test:
785     suffix: proj_tri_mdx_5P
786     requires: triangle
787     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
788     filter: grep -v marker | grep -v atomic | grep -v usage
789 
790   test:
791     suffix: proj_tri_3d
792     requires: ctetgen
793     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
794     filter: grep -v marker | grep -v atomic | grep -v usage
795 
796   test:
797     suffix: proj_tri_3d_2_faces
798     requires: ctetgen
799     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
800     filter: grep -v marker | grep -v atomic | grep -v usage
801 
802   test:
803     suffix: proj_tri_3d_5P
804     requires: ctetgen
805     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
806     filter: grep -v marker | grep -v atomic | grep -v usage
807 
808   test:
809     suffix: proj_tri_3d_mdx
810     requires: ctetgen
811     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
812     filter: grep -v marker | grep -v atomic | grep -v usage
813 
814   test:
815     suffix: proj_tri_3d_mdx_5P
816     requires: ctetgen
817     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
818     filter: grep -v marker | grep -v atomic | grep -v usage
819 
820   test:
821     suffix: proj_tri_3d_mdx_2_faces
822     requires: ctetgen
823     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
824     filter: grep -v marker | grep -v atomic | grep -v usage
825 
826   test:
827     suffix: proj_tri_3d_mdx_5P_2_faces
828     requires: ctetgen
829     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
830     filter: grep -v marker | grep -v atomic | grep -v usage
831 
832   test:
833     suffix: proj_quad_lsqr_scale
834     requires: !complex
835     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
836       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
837       -particlesPerCell 200 \
838       -ptof_pc_type lu  \
839       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none
840 
841   test:
842     suffix: proj_quad_lsqr_prec_scale
843     requires: !complex
844     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
845       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
846       -particlesPerCell 200 \
847       -ptof_pc_type lu  \
848       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec
849   test:
850     suffix: proj_shape_linear_tri_2d
851     requires: ctetgen triangle
852     args: -shape -dm_plex_box_faces 2,2\
853           -dm_view -sw_view\
854           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
855           -pc_type lu
856   test:
857     suffix: proj_shape_linear_quad_2d
858     requires: ctetgen triangle
859     args: -shape -dm_plex_simplex 0 -dm_plex_box_faces 2,2\
860           -dm_view -sw_view\
861           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
862           -pc_type lu
863   test:
864     suffix: proj_shape_linear_tri_3d
865     requires: ctetgen triangle
866     args: -shape -dm_plex_dim 3 -dm_plex_box_faces 2,2,2\
867           -dm_view -sw_view\
868           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
869           -pc_type lu
870   test:
871     suffix: proj_shape_linear_quad_3d
872     requires: ctetgen triangle
873     args: -shape -dm_plex_dim 3 -dm_plex_simplex 0\
874           -dm_plex_box_faces 2,2,2\
875           -dm_view -sw_view\
876           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
877           -pc_type lu
878 TEST*/
879