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