xref: /petsc/src/dm/impls/swarm/tests/ex2.c (revision fc3eb83c22be0eedfd551163d0b00bfc14ce0e8b)
180aecfd0Sjosephpu static char help[] = "Tests L2 projection with DMSwarm using delta function particles and deposition of linear shape.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscfe.h>
5c4762a1bSJed Brown #include <petscdmswarm.h>
6c4762a1bSJed Brown #include <petscds.h>
7c4762a1bSJed Brown #include <petscksp.h>
8c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h>
9c4762a1bSJed Brown typedef struct {
10832e2b15SMatthew G. Knepley   PetscReal L[3];             /* Dimensions of the mesh bounding box */
11c4762a1bSJed Brown   PetscInt  particlesPerCell; /* The number of partices per cell */
12c4762a1bSJed Brown   PetscReal particleRelDx;    /* Relative particle position perturbation compared to average cell diameter h */
13c4762a1bSJed Brown   PetscReal meshRelDx;        /* Relative vertex position perturbation compared to average cell diameter h */
14c4762a1bSJed Brown   PetscInt  k;                /* Mode number for test function */
15c4762a1bSJed Brown   PetscReal momentTol;        /* Tolerance for checking moment conservation */
1680aecfd0Sjosephpu   PetscBool shape;
17c4762a1bSJed Brown   PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
18c4762a1bSJed Brown } AppCtx;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */
21c4762a1bSJed Brown 
22d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
23d71ae5a4SJacob Faibussowitsch {
24c4762a1bSJed Brown   AppCtx  *ctx = (AppCtx *)a_ctx;
25c4762a1bSJed Brown   PetscInt d;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   u[0] = 0.0;
28832e2b15SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] += x[d] / (ctx->L[d]);
293ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32d71ae5a4SJacob Faibussowitsch static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
33d71ae5a4SJacob Faibussowitsch {
34c4762a1bSJed Brown   AppCtx  *ctx = (AppCtx *)a_ctx;
35c4762a1bSJed Brown   PetscInt d;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   u[0] = 1;
38832e2b15SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d]) * PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4);
393ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
42d71ae5a4SJacob Faibussowitsch static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown   AppCtx *ctx = (AppCtx *)a_ctx;
45c4762a1bSJed Brown 
46832e2b15SMatthew G. Knepley   u[0] = PetscSinScalar(2 * PETSC_PI * ctx->k * x[0] / (ctx->L[0]));
473ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   char      fstring[PETSC_MAX_PATH_LEN] = "linear";
53c4762a1bSJed Brown   PetscBool flag;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBeginUser;
56c4762a1bSJed Brown   options->particlesPerCell = 1;
57c4762a1bSJed Brown   options->k                = 1;
58c4762a1bSJed Brown   options->particleRelDx    = 1.e-20;
59c4762a1bSJed Brown   options->meshRelDx        = 1.e-20;
60c4762a1bSJed Brown   options->momentTol        = 100. * PETSC_MACHINE_EPSILON;
6180aecfd0Sjosephpu   options->shape            = PETSC_FALSE;
62c4762a1bSJed Brown 
63d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");
6480aecfd0Sjosephpu   PetscCall(PetscOptionsBool("-shape", "Test linear function projection", "ex2.c", options->shape, &options->shape, NULL));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL));
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL));
709566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(fstring, "linear", &flag));
71c4762a1bSJed Brown   if (flag) {
72c4762a1bSJed Brown     options->func = linear;
73c4762a1bSJed Brown   } else {
749566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(fstring, "sin", &flag));
75c4762a1bSJed Brown     if (flag) {
76c4762a1bSJed Brown       options->func = sinx;
77c4762a1bSJed Brown     } else {
789566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(fstring, "x2_x4", &flag));
79c4762a1bSJed Brown       options->func = x2_x4;
8028b400f6SJacob Faibussowitsch       PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown function %s", fstring);
81c4762a1bSJed Brown     }
82c4762a1bSJed Brown   }
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL));
84d0609cedSBarry Smith   PetscOptionsEnd();
85c4762a1bSJed Brown 
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89d71ae5a4SJacob Faibussowitsch static PetscErrorCode PerturbVertices(DM dm, AppCtx *user)
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown   PetscRandom  rnd;
92c4762a1bSJed Brown   PetscReal    interval = user->meshRelDx;
93c4762a1bSJed Brown   Vec          coordinates;
94c4762a1bSJed Brown   PetscScalar *coords;
95832e2b15SMatthew G. Knepley   PetscReal   *hh, low[3], high[3];
96832e2b15SMatthew G. Knepley   PetscInt     d, cdim, cEnd, N, p, bs;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   PetscFunctionBeginUser;
999566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dm, low, high));
1009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
1019566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
1029566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -interval, interval));
1039566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
1049566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1059566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
1069566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(cdim, &hh));
107832e2b15SMatthew G. Knepley   for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim);
1089566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(coordinates, &N));
1099566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(coordinates, &bs));
11063a3b9bcSJacob Faibussowitsch   PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim);
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
112c4762a1bSJed Brown   for (p = 0; p < N; p += cdim) {
113c4762a1bSJed Brown     PetscScalar *coord = &coords[p], value;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown     for (d = 0; d < cdim; ++d) {
1169566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rnd, &value));
117832e2b15SMatthew G. Knepley       coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d])));
118c4762a1bSJed Brown     }
119c4762a1bSJed Brown   }
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
1219566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1229566063dSJacob Faibussowitsch   PetscCall(PetscFree(hh));
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
127d71ae5a4SJacob Faibussowitsch {
128832e2b15SMatthew G. Knepley   PetscReal low[3], high[3];
129832e2b15SMatthew G. Knepley   PetscInt  cdim, d;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   PetscFunctionBeginUser;
1329566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1339566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1349566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1359566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
13630602db0SMatthew G. Knepley 
1379566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(*dm, &cdim));
1389566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(*dm, low, high));
139832e2b15SMatthew G. Knepley   for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d];
1409566063dSJacob Faibussowitsch   PetscCall(PerturbVertices(*dm, user));
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144d71ae5a4SJacob Faibussowitsch 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[])
145d71ae5a4SJacob Faibussowitsch {
146c4762a1bSJed Brown   g0[0] = 1.0;
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
150d71ae5a4SJacob Faibussowitsch {
151c4762a1bSJed Brown   PetscFE        fe;
152c4762a1bSJed Brown   PetscDS        ds;
153832e2b15SMatthew G. Knepley   DMPolytopeType ct;
154832e2b15SMatthew G. Knepley   PetscInt       dim, cStart;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBeginUser;
1579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
160*fc3eb83cSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PetscObjectComm((PetscObject)dm), dim, 1, ct, NULL, -1, &fe));
1619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
1629566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1639566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1649566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
165c4762a1bSJed Brown   /* Setup to form mass matrix */
1669566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1679566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL));
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
172d71ae5a4SJacob Faibussowitsch {
173c4762a1bSJed Brown   PetscRandom    rnd, rndp;
174c4762a1bSJed Brown   PetscReal      interval = user->particleRelDx;
175832e2b15SMatthew G. Knepley   DMPolytopeType ct;
176832e2b15SMatthew G. Knepley   PetscBool      simplex;
177c4762a1bSJed Brown   PetscScalar    value, *vals;
178c4762a1bSJed Brown   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
179c4762a1bSJed Brown   PetscInt      *cellid;
180832e2b15SMatthew G. Knepley   PetscInt       Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d;
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   PetscFunctionBeginUser;
1839566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1849566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
186832e2b15SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
1879566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1889566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1899566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
190c4762a1bSJed Brown 
1919566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
1929566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
1939566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
1949566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rndp));
1959566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rndp, -interval, interval));
1969566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rndp));
197c4762a1bSJed Brown 
1989566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
2019566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
2029566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell));
2039566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0));
2049566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
2059566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
2069566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
2079566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals));
208c4762a1bSJed Brown 
2099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
210c4762a1bSJed Brown   for (c = 0; c < Ncell; ++c) {
211c4762a1bSJed Brown     if (Np == 1) {
2129566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
213c4762a1bSJed Brown       cellid[c] = c;
214c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
215c4762a1bSJed Brown     } else {
2169566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
217c4762a1bSJed Brown       for (d = 0; d < dim; ++d) xi0[d] = -1.0;
218c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
219c4762a1bSJed Brown         const PetscInt n   = c * Np + p;
220c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
221c4762a1bSJed Brown 
222c4762a1bSJed Brown         cellid[n] = c;
2239371c9d4SSatish Balay         for (d = 0; d < dim; ++d) {
2249371c9d4SSatish Balay           PetscCall(PetscRandomGetValue(rnd, &value));
2259371c9d4SSatish Balay           refcoords[d] = PetscRealPart(value);
2269371c9d4SSatish Balay           sum += refcoords[d];
2279371c9d4SSatish Balay         }
2289371c9d4SSatish Balay         if (simplex && sum > 0.0)
2299371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
230c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
231c4762a1bSJed Brown       }
232c4762a1bSJed Brown     }
233c4762a1bSJed Brown   }
2349566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
235c4762a1bSJed Brown   for (c = 0; c < Ncell; ++c) {
236c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
237c4762a1bSJed Brown       const PetscInt n = c * Np + p;
238c4762a1bSJed Brown 
2399371c9d4SSatish Balay       for (d = 0; d < dim; ++d) {
2409371c9d4SSatish Balay         PetscCall(PetscRandomGetValue(rndp, &value));
2419371c9d4SSatish Balay         coords[n * dim + d] += PetscRealPart(value);
2429371c9d4SSatish Balay       }
2433ba16761SJacob Faibussowitsch       PetscCall(user->func(dim, 0.0, &coords[n * dim], 1, &vals[c], user));
244c4762a1bSJed Brown     }
245c4762a1bSJed Brown   }
2469566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
2479566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
2489566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals));
2499566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
2509566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rndp));
2519566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
2529566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
25680aecfd0Sjosephpu static PetscErrorCode CreateParticles_Shape(DM dm, DM *sw, AppCtx *user)
25780aecfd0Sjosephpu {
25880aecfd0Sjosephpu   PetscDS          prob;
25980aecfd0Sjosephpu   PetscFE          fe;
26080aecfd0Sjosephpu   PetscQuadrature  quad;
26180aecfd0Sjosephpu   PetscScalar     *vals;
26280aecfd0Sjosephpu   PetscReal       *v0, *J, *invJ, detJ, *coords, *xi0;
26380aecfd0Sjosephpu   PetscInt        *cellid;
26480aecfd0Sjosephpu   const PetscReal *qpoints, *qweights;
26580aecfd0Sjosephpu   PetscInt         cStart, cEnd, c, Nq, q, dim;
26680aecfd0Sjosephpu 
26780aecfd0Sjosephpu   PetscFunctionBeginUser;
26880aecfd0Sjosephpu   PetscCall(DMGetDimension(dm, &dim));
26980aecfd0Sjosephpu   PetscCall(DMGetDS(dm, &prob));
27080aecfd0Sjosephpu   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
27180aecfd0Sjosephpu   PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
27280aecfd0Sjosephpu   PetscCall(PetscFEGetQuadrature(fe, &quad));
27380aecfd0Sjosephpu   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, &qweights));
27480aecfd0Sjosephpu   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
27580aecfd0Sjosephpu   PetscCall(DMSetType(*sw, DMSWARM));
27680aecfd0Sjosephpu   PetscCall(DMSetDimension(*sw, dim));
27780aecfd0Sjosephpu   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
27880aecfd0Sjosephpu   PetscCall(DMSwarmSetCellDM(*sw, dm));
27980aecfd0Sjosephpu   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_REAL));
28080aecfd0Sjosephpu   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
28180aecfd0Sjosephpu   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Nq, 0));
28280aecfd0Sjosephpu   PetscCall(DMSetFromOptions(*sw));
28380aecfd0Sjosephpu   PetscCall(PetscMalloc4(dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
28480aecfd0Sjosephpu   for (c = 0; c < dim; c++) xi0[c] = -1.;
28580aecfd0Sjosephpu   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
28680aecfd0Sjosephpu   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
28780aecfd0Sjosephpu   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals));
28880aecfd0Sjosephpu   for (c = cStart; c < cEnd; ++c) {
28980aecfd0Sjosephpu     for (q = 0; q < Nq; ++q) {
29080aecfd0Sjosephpu       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
29180aecfd0Sjosephpu       cellid[c * Nq + q] = c;
29280aecfd0Sjosephpu       CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], &coords[(c * Nq + q) * dim]);
29380aecfd0Sjosephpu       PetscCall(linear(dim, 0.0, &coords[(c * Nq + q) * dim], 1, &vals[c * Nq + q], (void *)user));
29480aecfd0Sjosephpu       vals[c * Nq + q] *= qweights[q] * detJ;
29580aecfd0Sjosephpu     }
29680aecfd0Sjosephpu   }
29780aecfd0Sjosephpu   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
29880aecfd0Sjosephpu   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid));
29980aecfd0Sjosephpu   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals));
30080aecfd0Sjosephpu   PetscCall(PetscFree4(xi0, v0, J, invJ));
30180aecfd0Sjosephpu   PetscCall(DMSwarmMigrate(*sw, PETSC_FALSE));
30280aecfd0Sjosephpu   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
30380aecfd0Sjosephpu   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
30480aecfd0Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
30580aecfd0Sjosephpu }
30680aecfd0Sjosephpu 
307d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user)
308d71ae5a4SJacob Faibussowitsch {
309c4762a1bSJed Brown   DM                 dm;
310c4762a1bSJed Brown   const PetscReal   *coords;
311c4762a1bSJed Brown   const PetscScalar *w;
312c4762a1bSJed Brown   PetscReal          mom[3] = {0.0, 0.0, 0.0};
313c4762a1bSJed Brown   PetscInt           cell, cStart, cEnd, dim;
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   PetscFunctionBeginUser;
3169566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
3179566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dm));
3189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3199566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
3209566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
3219566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w));
322c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
323c4762a1bSJed Brown     PetscInt *pidx;
324c4762a1bSJed Brown     PetscInt  Np, p, d;
325c4762a1bSJed Brown 
3269566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
327c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
328c4762a1bSJed Brown       const PetscInt   idx = pidx[p];
329c4762a1bSJed Brown       const PetscReal *c   = &coords[idx * dim];
330c4762a1bSJed Brown 
331c4762a1bSJed Brown       mom[0] += PetscRealPart(w[idx]);
332c4762a1bSJed Brown       mom[1] += PetscRealPart(w[idx]) * c[0];
333c4762a1bSJed Brown       for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d];
334c4762a1bSJed Brown     }
3359566063dSJacob Faibussowitsch     PetscCall(PetscFree(pidx));
336c4762a1bSJed Brown   }
3379566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
3389566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
3399566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
340712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown 
344d71ae5a4SJacob Faibussowitsch 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[])
345d71ae5a4SJacob Faibussowitsch {
346c4762a1bSJed Brown   f0[0] = u[0];
347c4762a1bSJed Brown }
348c4762a1bSJed Brown 
349d71ae5a4SJacob Faibussowitsch 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[])
350d71ae5a4SJacob Faibussowitsch {
351c4762a1bSJed Brown   f0[0] = x[0] * u[0];
352c4762a1bSJed Brown }
353c4762a1bSJed Brown 
354d71ae5a4SJacob Faibussowitsch 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[])
355d71ae5a4SJacob Faibussowitsch {
356c4762a1bSJed Brown   PetscInt d;
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   f0[0] = 0.0;
359c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0];
360c4762a1bSJed Brown }
361c4762a1bSJed Brown 
362d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
363d71ae5a4SJacob Faibussowitsch {
364c4762a1bSJed Brown   PetscDS     prob;
365c4762a1bSJed Brown   PetscScalar mom;
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   PetscFunctionBeginUser;
3689566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
3699566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_1));
3709566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
371c4762a1bSJed Brown   moments[0] = PetscRealPart(mom);
3729566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_x));
3739566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
374c4762a1bSJed Brown   moments[1] = PetscRealPart(mom);
3759566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_r2));
3769566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
377c4762a1bSJed Brown   moments[2] = PetscRealPart(mom);
3783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
379c4762a1bSJed Brown }
380c4762a1bSJed Brown 
381*fc3eb83cSMatthew G. Knepley static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, Vec fhat, AppCtx *user)
382d71ae5a4SJacob Faibussowitsch {
383*fc3eb83cSMatthew G. Knepley   const char *fieldnames[1] = {"w_q"};
384*fc3eb83cSMatthew G. Knepley   Vec         fields[1]     = {fhat};
385*fc3eb83cSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
386*fc3eb83cSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   PetscFunctionBeginUser;
389*fc3eb83cSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, 1, fieldnames, fields, SCATTER_FORWARD));
390c4762a1bSJed Brown 
391c4762a1bSJed Brown   /* Check moments of field */
3929566063dSJacob Faibussowitsch   PetscCall(computeParticleMoments(sw, pmoments, user));
3939566063dSJacob Faibussowitsch   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
394*fc3eb83cSMatthew G. Knepley   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
395*fc3eb83cSMatthew G. Knepley   for (PetscInt m = 0; m < 3; ++m) {
396*fc3eb83cSMatthew G. Knepley     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]),
397*fc3eb83cSMatthew G. Knepley                user->momentTol);
398c4762a1bSJed Brown   }
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
400c4762a1bSJed Brown }
401c4762a1bSJed Brown 
402*fc3eb83cSMatthew G. Knepley static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, Vec fhat, AppCtx *user)
403d71ae5a4SJacob Faibussowitsch {
404*fc3eb83cSMatthew G. Knepley   const char *fieldnames[1] = {"w_q"};
405*fc3eb83cSMatthew G. Knepley   Vec         fields[1]     = {fhat};
406*fc3eb83cSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
407*fc3eb83cSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   PetscFunctionBeginUser;
410*fc3eb83cSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, 1, fieldnames, fields, SCATTER_REVERSE));
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   /* Check moments */
4139566063dSJacob Faibussowitsch   PetscCall(computeParticleMoments(sw, pmoments, user));
4149566063dSJacob Faibussowitsch   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
415*fc3eb83cSMatthew G. Knepley   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
416*fc3eb83cSMatthew G. Knepley   for (PetscInt m = 0; m < 3; ++m) {
417*fc3eb83cSMatthew G. Knepley     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]),
418*fc3eb83cSMatthew G. Knepley                user->momentTol);
419c4762a1bSJed Brown   }
4203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
421c4762a1bSJed Brown }
422c4762a1bSJed Brown 
423c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */
424d71ae5a4SJacob Faibussowitsch static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC)
425d71ae5a4SJacob Faibussowitsch {
426c4762a1bSJed Brown   DM_Plex         *mesh  = (DM_Plex *)dm->data;
427c4762a1bSJed Brown   PetscInt         debug = mesh->printFEM;
428c4762a1bSJed Brown   DM               dmC;
429c4762a1bSJed Brown   PetscSection     section;
430c4762a1bSJed Brown   PetscQuadrature  quad = NULL;
431c4762a1bSJed Brown   PetscScalar     *interpolant, *gradsum;
432c4762a1bSJed Brown   PetscFEGeom      fegeom;
433c4762a1bSJed Brown   PetscReal       *coords;
434c4762a1bSJed Brown   const PetscReal *quadPoints, *quadWeights;
435c4762a1bSJed Brown   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   PetscFunctionBegin;
4389566063dSJacob Faibussowitsch   PetscCall(VecGetDM(locC, &dmC));
4399566063dSJacob Faibussowitsch   PetscCall(VecSet(locC, 0.0));
4409566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
442c4762a1bSJed Brown   fegeom.dimEmbed = coordDim;
4439566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
4449566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &numFields));
445c4762a1bSJed Brown   for (field = 0; field < numFields; ++field) {
446c4762a1bSJed Brown     PetscObject  obj;
447c4762a1bSJed Brown     PetscClassId id;
448c4762a1bSJed Brown     PetscInt     Nc;
449c4762a1bSJed Brown 
4509566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, field, NULL, &obj));
4519566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
452c4762a1bSJed Brown     if (id == PETSCFE_CLASSID) {
453c4762a1bSJed Brown       PetscFE fe = (PetscFE)obj;
454c4762a1bSJed Brown 
4559566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature(fe, &quad));
4569566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(fe, &Nc));
457c4762a1bSJed Brown     } else if (id == PETSCFV_CLASSID) {
458c4762a1bSJed Brown       PetscFV fv = (PetscFV)obj;
459c4762a1bSJed Brown 
4609566063dSJacob Faibussowitsch       PetscCall(PetscFVGetQuadrature(fv, &quad));
4619566063dSJacob Faibussowitsch       PetscCall(PetscFVGetNumComponents(fv, &Nc));
46263a3b9bcSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
463c4762a1bSJed Brown     numComponents += Nc;
464c4762a1bSJed Brown   }
4659566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
46663a3b9bcSJacob Faibussowitsch   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
4679566063dSJacob Faibussowitsch   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));
4689566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4699566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
470c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
471c4762a1bSJed Brown     PetscInt *star = NULL;
472c4762a1bSJed Brown     PetscInt  starSize, st, d, fc;
473c4762a1bSJed Brown 
4749566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(gradsum, coordDim * numComponents));
4759566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
476c4762a1bSJed Brown     for (st = 0; st < starSize * 2; st += 2) {
477c4762a1bSJed Brown       const PetscInt cell = star[st];
478c4762a1bSJed Brown       PetscScalar   *grad = &gradsum[coordDim * numComponents];
479c4762a1bSJed Brown       PetscScalar   *x    = NULL;
480c4762a1bSJed Brown 
481c4762a1bSJed Brown       if ((cell < cStart) || (cell >= cEnd)) continue;
4829566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
4839566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
484c4762a1bSJed Brown       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
485c4762a1bSJed Brown         PetscObject  obj;
486c4762a1bSJed Brown         PetscClassId id;
487c4762a1bSJed Brown         PetscInt     Nb, Nc, q, qc = 0;
488c4762a1bSJed Brown 
4899566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(grad, coordDim * numComponents));
4909566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, field, NULL, &obj));
4919566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetClassId(obj, &id));
4929371c9d4SSatish Balay         if (id == PETSCFE_CLASSID) {
4939371c9d4SSatish Balay           PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
4949371c9d4SSatish Balay           PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
4959371c9d4SSatish Balay         } else if (id == PETSCFV_CLASSID) {
4969371c9d4SSatish Balay           PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
4979371c9d4SSatish Balay           Nb = 1;
4989371c9d4SSatish Balay         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
499c4762a1bSJed Brown         for (q = 0; q < Nq; ++q) {
50063a3b9bcSJacob Faibussowitsch           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);
5019566063dSJacob Faibussowitsch           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant));
50263a3b9bcSJacob Faibussowitsch           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
503c4762a1bSJed Brown           for (fc = 0; fc < Nc; ++fc) {
504c4762a1bSJed Brown             const PetscReal wt = quadWeights[q * qNc + qc + fc];
505c4762a1bSJed Brown 
506c4762a1bSJed Brown             for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q];
507c4762a1bSJed Brown           }
508c4762a1bSJed Brown         }
509c4762a1bSJed Brown         fieldOffset += Nb;
510c4762a1bSJed Brown       }
5119566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
512c4762a1bSJed Brown       for (fc = 0; fc < numComponents; ++fc) {
513ad540459SPierre Jolivet         for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d];
514c4762a1bSJed Brown       }
515c4762a1bSJed Brown       if (debug) {
51663a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell));
517c4762a1bSJed Brown         for (fc = 0; fc < numComponents; ++fc) {
518c4762a1bSJed Brown           for (d = 0; d < coordDim; ++d) {
5199566063dSJacob Faibussowitsch             if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
5209566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d])));
521c4762a1bSJed Brown           }
522c4762a1bSJed Brown         }
5239566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
524c4762a1bSJed Brown       }
525c4762a1bSJed Brown     }
5269566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
5279566063dSJacob Faibussowitsch     PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
528c4762a1bSJed Brown   }
5299566063dSJacob Faibussowitsch   PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531c4762a1bSJed Brown }
532c4762a1bSJed Brown 
533d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user)
534d71ae5a4SJacob Faibussowitsch {
535c4762a1bSJed Brown   MPI_Comm  comm;
536c4762a1bSJed Brown   KSP       ksp;
537c4762a1bSJed Brown   Mat       M;                  /* FEM mass matrix */
538c4762a1bSJed Brown   Mat       M_p;                /* Particle mass matrix */
539c4762a1bSJed Brown   Vec       f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */
540c4762a1bSJed Brown   PetscReal pmoments[3];        /* \int f, \int x f, \int r^2 f */
541c4762a1bSJed Brown   PetscReal fmoments[3];        /* \int \hat f, \int x \hat f, \int r^2 \hat f */
542c4762a1bSJed Brown   PetscInt  m;
543c4762a1bSJed Brown 
544c4762a1bSJed Brown   PetscFunctionBeginUser;
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5469566063dSJacob Faibussowitsch   PetscCall(KSPCreate(comm, &ksp));
5479566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(ksp, "ptof_"));
5489566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &fhat));
5499566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &rhs));
5509566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &grad));
551c4762a1bSJed Brown 
5529566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
5539566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
554c4762a1bSJed Brown 
555c4762a1bSJed Brown   /* make particle weight vector */
5569566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
557c4762a1bSJed Brown 
558c4762a1bSJed Brown   /* create matrix RHS vector */
5599566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(M_p, f, rhs));
5609566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
5619566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs"));
5629566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
563c4762a1bSJed Brown 
5649566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &M));
5659566063dSJacob Faibussowitsch   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
566c4762a1bSJed Brown 
5679566063dSJacob Faibussowitsch   PetscCall(InterpolateGradient(dm, fhat, grad));
568c4762a1bSJed Brown 
5699566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
5709566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, M, M));
5719566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
5729566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, rhs, grad));
5739566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat"));
5749566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
575c4762a1bSJed Brown 
576c4762a1bSJed Brown   /* Check moments of field */
5779566063dSJacob Faibussowitsch   PetscCall(computeParticleMoments(sw, pmoments, user));
5789566063dSJacob Faibussowitsch   PetscCall(computeFEMMoments(dm, grad, fmoments, user));
5799566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
580c4762a1bSJed Brown   for (m = 0; m < 3; ++m) {
5811dca8a05SBarry Smith     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);
582c4762a1bSJed Brown   }
583c4762a1bSJed Brown 
5849566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
5859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M));
5869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M_p));
5879566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &fhat));
5889566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &rhs));
5899566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &grad));
590c4762a1bSJed Brown 
5913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
592c4762a1bSJed Brown }
593c4762a1bSJed Brown 
59480aecfd0Sjosephpu static PetscErrorCode TestL2Projection_Shape(DM dm, DM sw, AppCtx *user)
59580aecfd0Sjosephpu {
59680aecfd0Sjosephpu   PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
59780aecfd0Sjosephpu   KSP       ksp;
59880aecfd0Sjosephpu   Mat       mass;
59980aecfd0Sjosephpu   Vec       u, rhs, uproj;
60080aecfd0Sjosephpu   void     *ctxs[1];
60180aecfd0Sjosephpu   PetscReal error;
60280aecfd0Sjosephpu 
60380aecfd0Sjosephpu   PetscFunctionBeginUser;
60480aecfd0Sjosephpu   ctxs[0]  = (void *)user;
60580aecfd0Sjosephpu   funcs[0] = linear;
60680aecfd0Sjosephpu 
60780aecfd0Sjosephpu   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &u));
60880aecfd0Sjosephpu   PetscCall(VecViewFromOptions(u, NULL, "-f_view"));
60980aecfd0Sjosephpu   PetscCall(DMGetGlobalVector(dm, &rhs));
61080aecfd0Sjosephpu   PetscCall(DMCreateMassMatrix(sw, dm, &mass));
61180aecfd0Sjosephpu   PetscCall(MatMultTranspose(mass, u, rhs));
61280aecfd0Sjosephpu   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
61380aecfd0Sjosephpu   PetscCall(MatDestroy(&mass));
61480aecfd0Sjosephpu   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &u));
61580aecfd0Sjosephpu   PetscCall(DMGetGlobalVector(dm, &uproj));
61680aecfd0Sjosephpu   PetscCall(DMCreateMatrix(dm, &mass));
61780aecfd0Sjosephpu   PetscCall(DMPlexSNESComputeJacobianFEM(dm, uproj, mass, mass, user));
61880aecfd0Sjosephpu   PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view"));
61980aecfd0Sjosephpu   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
62080aecfd0Sjosephpu   PetscCall(KSPSetOperators(ksp, mass, mass));
62180aecfd0Sjosephpu   PetscCall(KSPSetFromOptions(ksp));
62280aecfd0Sjosephpu   PetscCall(KSPSolve(ksp, rhs, uproj));
62380aecfd0Sjosephpu   PetscCall(KSPDestroy(&ksp));
62480aecfd0Sjosephpu   PetscCall(PetscObjectSetName((PetscObject)uproj, "Full Projection"));
62580aecfd0Sjosephpu   PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view"));
62680aecfd0Sjosephpu   PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, uproj, &error));
62780aecfd0Sjosephpu   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error));
62880aecfd0Sjosephpu   PetscCall(MatDestroy(&mass));
62980aecfd0Sjosephpu   PetscCall(DMRestoreGlobalVector(dm, &rhs));
63080aecfd0Sjosephpu   PetscCall(DMRestoreGlobalVector(dm, &uproj));
63180aecfd0Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
63280aecfd0Sjosephpu }
63380aecfd0Sjosephpu 
634d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
635d71ae5a4SJacob Faibussowitsch {
636c4762a1bSJed Brown   MPI_Comm comm;
637c4762a1bSJed Brown   DM       dm, sw;
638c4762a1bSJed Brown   AppCtx   user;
639c4762a1bSJed Brown 
640327415f7SBarry Smith   PetscFunctionBeginUser;
6419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
642c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
6439566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
6449566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
6459566063dSJacob Faibussowitsch   PetscCall(CreateFEM(dm, &user));
64680aecfd0Sjosephpu   if (!user.shape) {
647*fc3eb83cSMatthew G. Knepley     Vec fhat;
648*fc3eb83cSMatthew G. Knepley 
6499566063dSJacob Faibussowitsch     PetscCall(CreateParticles(dm, &sw, &user));
650*fc3eb83cSMatthew G. Knepley     PetscCall(DMGetGlobalVector(dm, &fhat));
651*fc3eb83cSMatthew G. Knepley     PetscCall(TestL2ProjectionParticlesToField(dm, sw, fhat, &user));
652*fc3eb83cSMatthew G. Knepley     PetscCall(TestL2ProjectionFieldToParticles(dm, sw, fhat, &user));
6539566063dSJacob Faibussowitsch     PetscCall(TestFieldGradientProjection(dm, sw, &user));
654*fc3eb83cSMatthew G. Knepley     PetscCall(DMRestoreGlobalVector(dm, &fhat));
65580aecfd0Sjosephpu   } else {
65680aecfd0Sjosephpu     PetscCall(CreateParticles_Shape(dm, &sw, &user));
65780aecfd0Sjosephpu     PetscCall(TestL2Projection_Shape(dm, sw, &user));
65880aecfd0Sjosephpu   }
6599566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
6609566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
6619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
662b122ec5aSJacob Faibussowitsch   return 0;
663c4762a1bSJed Brown }
664c4762a1bSJed Brown 
665c4762a1bSJed Brown /*TEST
666c4762a1bSJed Brown 
667832e2b15SMatthew G. Knepley   # Swarm does not handle complex or quad
668832e2b15SMatthew G. Knepley   build:
669832e2b15SMatthew G. Knepley     requires: !complex double
670832e2b15SMatthew G. Knepley 
671c4762a1bSJed Brown   test:
672c4762a1bSJed Brown     suffix: proj_tri_0
673832e2b15SMatthew G. Knepley     requires: triangle
674832e2b15SMatthew G. Knepley     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
675c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
676c4762a1bSJed Brown 
677c4762a1bSJed Brown   test:
678c4762a1bSJed Brown     suffix: proj_tri_2_faces
679832e2b15SMatthew G. Knepley     requires: triangle
680832e2b15SMatthew G. Knepley     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
681c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
682c4762a1bSJed Brown 
683c4762a1bSJed Brown   test:
684c4762a1bSJed Brown     suffix: proj_quad_0
685832e2b15SMatthew G. Knepley     requires: triangle
68630602db0SMatthew G. Knepley     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
687c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
688c4762a1bSJed Brown 
689c4762a1bSJed Brown   test:
690c4762a1bSJed Brown     suffix: proj_quad_2_faces
691832e2b15SMatthew G. Knepley     requires: triangle
69230602db0SMatthew G. Knepley     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
693c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
694c4762a1bSJed Brown 
695c4762a1bSJed Brown   test:
696c4762a1bSJed Brown     suffix: proj_tri_5P
697832e2b15SMatthew G. Knepley     requires: triangle
698832e2b15SMatthew G. Knepley     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
699c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
700c4762a1bSJed Brown 
701c4762a1bSJed Brown   test:
702c4762a1bSJed Brown     suffix: proj_quad_5P
703832e2b15SMatthew G. Knepley     requires: triangle
70430602db0SMatthew G. Knepley     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
705c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
706c4762a1bSJed Brown 
707c4762a1bSJed Brown   test:
708c4762a1bSJed Brown     suffix: proj_tri_mdx
709832e2b15SMatthew G. Knepley     requires: triangle
710832e2b15SMatthew G. Knepley     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
711c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
712c4762a1bSJed Brown 
713c4762a1bSJed Brown   test:
714c4762a1bSJed Brown     suffix: proj_tri_mdx_5P
715832e2b15SMatthew G. Knepley     requires: triangle
716832e2b15SMatthew G. Knepley     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
717c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
718c4762a1bSJed Brown 
719c4762a1bSJed Brown   test:
720c4762a1bSJed Brown     suffix: proj_tri_3d
721832e2b15SMatthew G. Knepley     requires: ctetgen
72230602db0SMatthew G. Knepley     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
723c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
724c4762a1bSJed Brown 
725c4762a1bSJed Brown   test:
726c4762a1bSJed Brown     suffix: proj_tri_3d_2_faces
727832e2b15SMatthew G. Knepley     requires: ctetgen
72830602db0SMatthew G. Knepley     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
729c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
730c4762a1bSJed Brown 
731c4762a1bSJed Brown   test:
732c4762a1bSJed Brown     suffix: proj_tri_3d_5P
733832e2b15SMatthew G. Knepley     requires: ctetgen
73430602db0SMatthew G. Knepley     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
735c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
736c4762a1bSJed Brown 
737c4762a1bSJed Brown   test:
738c4762a1bSJed Brown     suffix: proj_tri_3d_mdx
739832e2b15SMatthew G. Knepley     requires: ctetgen
74030602db0SMatthew G. Knepley     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
741c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
742c4762a1bSJed Brown 
743c4762a1bSJed Brown   test:
744c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_5P
745832e2b15SMatthew G. Knepley     requires: ctetgen
74630602db0SMatthew G. Knepley     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
747c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
748c4762a1bSJed Brown 
749c4762a1bSJed Brown   test:
750c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_2_faces
751832e2b15SMatthew G. Knepley     requires: ctetgen
75230602db0SMatthew G. Knepley     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
753c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
754c4762a1bSJed Brown 
755c4762a1bSJed Brown   test:
756c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_5P_2_faces
757832e2b15SMatthew G. Knepley     requires: ctetgen
75830602db0SMatthew G. Knepley     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
759c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
760c4762a1bSJed Brown 
761832e2b15SMatthew G. Knepley   test:
762832e2b15SMatthew G. Knepley     suffix: proj_quad_lsqr_scale
763832e2b15SMatthew G. Knepley     requires: !complex
76430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
765832e2b15SMatthew G. Knepley       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
766832e2b15SMatthew G. Knepley       -particlesPerCell 200 \
767832e2b15SMatthew G. Knepley       -ptof_pc_type lu  \
768832e2b15SMatthew G. Knepley       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none
769832e2b15SMatthew G. Knepley 
770832e2b15SMatthew G. Knepley   test:
771832e2b15SMatthew G. Knepley     suffix: proj_quad_lsqr_prec_scale
772832e2b15SMatthew G. Knepley     requires: !complex
77330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
774832e2b15SMatthew G. Knepley       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
775832e2b15SMatthew G. Knepley       -particlesPerCell 200 \
776832e2b15SMatthew G. Knepley       -ptof_pc_type lu  \
777*fc3eb83cSMatthew G. Knepley       -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
77880aecfd0Sjosephpu   test:
77980aecfd0Sjosephpu     suffix: proj_shape_linear_tri_2d
78080aecfd0Sjosephpu     requires: ctetgen triangle
78180aecfd0Sjosephpu     args: -shape -dm_plex_box_faces 2,2\
78280aecfd0Sjosephpu           -dm_view -sw_view\
78380aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
78480aecfd0Sjosephpu           -pc_type lu
78580aecfd0Sjosephpu   test:
78680aecfd0Sjosephpu     suffix: proj_shape_linear_quad_2d
78780aecfd0Sjosephpu     requires: ctetgen triangle
78880aecfd0Sjosephpu     args: -shape -dm_plex_simplex 0 -dm_plex_box_faces 2,2\
78980aecfd0Sjosephpu           -dm_view -sw_view\
79080aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
79180aecfd0Sjosephpu           -pc_type lu
79280aecfd0Sjosephpu   test:
79380aecfd0Sjosephpu     suffix: proj_shape_linear_tri_3d
79480aecfd0Sjosephpu     requires: ctetgen triangle
79580aecfd0Sjosephpu     args: -shape -dm_plex_dim 3 -dm_plex_box_faces 2,2,2\
79680aecfd0Sjosephpu           -dm_view -sw_view\
79780aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
79880aecfd0Sjosephpu           -pc_type lu
79980aecfd0Sjosephpu   test:
80080aecfd0Sjosephpu     suffix: proj_shape_linear_quad_3d
80180aecfd0Sjosephpu     requires: ctetgen triangle
80280aecfd0Sjosephpu     args: -shape -dm_plex_dim 3 -dm_plex_simplex 0\
80380aecfd0Sjosephpu           -dm_plex_box_faces 2,2,2\
80480aecfd0Sjosephpu           -dm_view -sw_view\
80580aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
80680aecfd0Sjosephpu           -pc_type lu
807c4762a1bSJed Brown TEST*/
808