xref: /petsc/src/dm/impls/swarm/tests/ex2.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Tests L2 projection with DMSwarm using delta function particles.\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 */
16832e2b15SMatthew G. Knepley   PetscBool useBlockDiagPrec;                 /* Use the block diagonal of the normal equations as a preconditioner */
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 
22c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
23c4762a1bSJed Brown {
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]);
29c4762a1bSJed Brown   return 0;
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32c4762a1bSJed Brown static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
33c4762a1bSJed Brown {
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);
39c4762a1bSJed Brown   return 0;
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
42c4762a1bSJed Brown static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx)
43c4762a1bSJed Brown {
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]));
47c4762a1bSJed Brown   return 0;
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   char           fstring[PETSC_MAX_PATH_LEN] = "linear";
53c4762a1bSJed Brown   PetscBool      flag;
54c4762a1bSJed Brown   PetscErrorCode ierr;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   PetscFunctionBeginUser;
57c4762a1bSJed Brown   options->particlesPerCell = 1;
58c4762a1bSJed Brown   options->k                = 1;
59c4762a1bSJed Brown   options->particleRelDx    = 1.e-20;
60c4762a1bSJed Brown   options->meshRelDx        = 1.e-20;
61c4762a1bSJed Brown   options->momentTol        = 100.*PETSC_MACHINE_EPSILON;
62832e2b15SMatthew G. Knepley   options->useBlockDiagPrec = PETSC_FALSE;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");CHKERRQ(ierr);
655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(fstring, "linear", &flag));
71c4762a1bSJed Brown   if (flag) {
72c4762a1bSJed Brown     options->func = linear;
73c4762a1bSJed Brown   } else {
745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(fstring, "sin", &flag));
75c4762a1bSJed Brown     if (flag) {
76c4762a1bSJed Brown       options->func = sinx;
77c4762a1bSJed Brown     } else {
785f80ce2aSJacob Faibussowitsch       CHKERRQ(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   }
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-block_diag_prec", "Use the block diagonal of the normal equations to precondition the particle projection", "ex2.c", options->useBlockDiagPrec, &options->useBlockDiagPrec, NULL));
851e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   PetscFunctionReturn(0);
88c4762a1bSJed Brown }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown static PetscErrorCode PerturbVertices(DM dm, AppCtx *user)
91c4762a1bSJed Brown {
92c4762a1bSJed Brown   PetscRandom    rnd;
93c4762a1bSJed Brown   PetscReal      interval = user->meshRelDx;
94c4762a1bSJed Brown   Vec            coordinates;
95c4762a1bSJed Brown   PetscScalar   *coords;
96832e2b15SMatthew G. Knepley   PetscReal     *hh, low[3], high[3];
97832e2b15SMatthew G. Knepley   PetscInt       d, cdim, cEnd, N, p, bs;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   PetscFunctionBeginUser;
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetBoundingBox(dm, low, high));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rnd, -interval, interval));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rnd));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dm, &cdim));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(cdim,&hh));
108832e2b15SMatthew G. Knepley   for (d = 0; d < cdim; ++d) hh[d] = (user->L[d])/PetscPowReal(cEnd, 1./cdim);
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(coordinates, &N));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetBlockSize(coordinates, &bs));
1112c71b3e2SJacob Faibussowitsch   PetscCheckFalse(bs != cdim,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %D != %D", bs, cdim);
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(coordinates, &coords));
113c4762a1bSJed Brown   for (p = 0; p < N; p += cdim) {
114c4762a1bSJed Brown     PetscScalar *coord = &coords[p], value;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     for (d = 0; d < cdim; ++d) {
1175f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rnd, &value));
118832e2b15SMatthew G. Knepley       coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value*hh[d])));
119c4762a1bSJed Brown     }
120c4762a1bSJed Brown   }
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(coordinates, &coords));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(hh));
124c4762a1bSJed Brown   PetscFunctionReturn(0);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
128c4762a1bSJed Brown {
129832e2b15SMatthew G. Knepley   PetscReal      low[3], high[3];
130832e2b15SMatthew G. Knepley   PetscInt       cdim, d;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   PetscFunctionBeginUser;
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
13730602db0SMatthew G. Knepley 
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(*dm, &cdim));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetBoundingBox(*dm, low, high));
140832e2b15SMatthew G. Knepley   for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d];
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(PerturbVertices(*dm, user));
142c4762a1bSJed Brown   PetscFunctionReturn(0);
143c4762a1bSJed Brown }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux,
146c4762a1bSJed Brown                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
147c4762a1bSJed Brown                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
148c4762a1bSJed Brown                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
149c4762a1bSJed Brown {
150c4762a1bSJed Brown   g0[0] = 1.0;
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
154c4762a1bSJed Brown {
155c4762a1bSJed Brown   PetscFE        fe;
156c4762a1bSJed Brown   PetscDS        ds;
157832e2b15SMatthew G. Knepley   DMPolytopeType ct;
158832e2b15SMatthew G. Knepley   PetscBool      simplex;
159832e2b15SMatthew G. Knepley   PetscInt       dim, cStart;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   PetscFunctionBeginUser;
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(dm, cStart, &ct));
165832e2b15SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, NULL, -1, &fe));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, "fe"));
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
171c4762a1bSJed Brown   /* Setup to form mass matrix */
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL));
174c4762a1bSJed Brown   PetscFunctionReturn(0);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
178c4762a1bSJed Brown {
179c4762a1bSJed Brown   PetscRandom    rnd, rndp;
180c4762a1bSJed Brown   PetscReal      interval = user->particleRelDx;
181832e2b15SMatthew G. Knepley   DMPolytopeType ct;
182832e2b15SMatthew G. Knepley   PetscBool      simplex;
183c4762a1bSJed Brown   PetscScalar    value, *vals;
184c4762a1bSJed Brown   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
185c4762a1bSJed Brown   PetscInt      *cellid;
186832e2b15SMatthew G. Knepley   PetscInt       Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   PetscFunctionBeginUser;
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(dm, cStart, &ct));
192832e2b15SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*sw, DMSWARM));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(*sw, dim));
196c4762a1bSJed Brown 
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rnd));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rndp));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rndp, -interval, interval));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rndp));
203c4762a1bSJed Brown 
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(*sw, dm));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(*sw));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*sw));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **) &vals));
214c4762a1bSJed Brown 
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
216c4762a1bSJed Brown   for (c = 0; c < Ncell; ++c) {
217c4762a1bSJed Brown     if (Np == 1) {
2185f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
219c4762a1bSJed Brown       cellid[c] = c;
220c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
221c4762a1bSJed Brown     } else {
2225f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
223c4762a1bSJed Brown       for (d = 0; d < dim; ++d) xi0[d] = -1.0;
224c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
225c4762a1bSJed Brown         const PetscInt n   = c*Np + p;
226c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
227c4762a1bSJed Brown 
228c4762a1bSJed Brown         cellid[n] = c;
2295f80ce2aSJacob Faibussowitsch         for (d = 0; d < dim; ++d) {CHKERRQ(PetscRandomGetValue(rnd, &value)); refcoords[d] = PetscRealPart(value); sum += refcoords[d];}
230832e2b15SMatthew G. Knepley         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
231c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
232c4762a1bSJed Brown       }
233c4762a1bSJed Brown     }
234c4762a1bSJed Brown   }
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(centroid, xi0, v0, J, invJ));
236c4762a1bSJed Brown   for (c = 0; c < Ncell; ++c) {
237c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
238c4762a1bSJed Brown       const PetscInt n = c*Np + p;
239c4762a1bSJed Brown 
2405f80ce2aSJacob Faibussowitsch       for (d = 0; d < dim; ++d) {CHKERRQ(PetscRandomGetValue(rndp, &value)); coords[n*dim+d] += PetscRealPart(value);}
241c4762a1bSJed Brown       user->func(dim, 0.0, &coords[n*dim], 1, &vals[c], user);
242c4762a1bSJed Brown     }
243c4762a1bSJed Brown   }
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **) &vals));
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rndp));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles"));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*sw, NULL, "-sw_view"));
251c4762a1bSJed Brown   PetscFunctionReturn(0);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
254c4762a1bSJed Brown static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user)
255c4762a1bSJed Brown {
256c4762a1bSJed Brown   DM                 dm;
257c4762a1bSJed Brown   const PetscReal   *coords;
258c4762a1bSJed Brown   const PetscScalar *w;
259c4762a1bSJed Brown   PetscReal          mom[3] = {0.0, 0.0, 0.0};
260c4762a1bSJed Brown   PetscInt           cell, cStart, cEnd, dim;
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   PetscFunctionBeginUser;
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(sw, &dim));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(sw, &dm));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortGetAccess(sw));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **) &w));
269c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
270c4762a1bSJed Brown     PetscInt *pidx;
271c4762a1bSJed Brown     PetscInt  Np, p, d;
272c4762a1bSJed Brown 
2735f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
274c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
275c4762a1bSJed Brown       const PetscInt   idx = pidx[p];
276c4762a1bSJed Brown       const PetscReal *c   = &coords[idx*dim];
277c4762a1bSJed Brown 
278c4762a1bSJed Brown       mom[0] += PetscRealPart(w[idx]);
279c4762a1bSJed Brown       mom[1] += PetscRealPart(w[idx]) * c[0];
280c4762a1bSJed Brown       for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d]*c[d];
281c4762a1bSJed Brown     }
2825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(pidx));
283c4762a1bSJed Brown   }
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **) &w));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSortRestoreAccess(sw));
2875f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) sw)));
288c4762a1bSJed Brown   PetscFunctionReturn(0);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
291c4762a1bSJed Brown static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
292c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
293c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
294c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
295c4762a1bSJed Brown {
296c4762a1bSJed Brown   f0[0] = u[0];
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux,
300c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
301c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
302c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
303c4762a1bSJed Brown {
304c4762a1bSJed Brown   f0[0] = x[0]*u[0];
305c4762a1bSJed Brown }
306c4762a1bSJed Brown 
307c4762a1bSJed Brown static void f0_r2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
308c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
309c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
310c4762a1bSJed Brown                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
311c4762a1bSJed Brown {
312c4762a1bSJed Brown   PetscInt d;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   f0[0] = 0.0;
315c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d])*u[0];
316c4762a1bSJed Brown }
317c4762a1bSJed Brown 
318c4762a1bSJed Brown static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
319c4762a1bSJed Brown {
320c4762a1bSJed Brown   PetscDS        prob;
321c4762a1bSJed Brown   PetscScalar    mom;
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   PetscFunctionBeginUser;
3245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &prob));
3255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetObjective(prob, 0, &f0_1));
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexComputeIntegralFEM(dm, u, &mom, user));
327c4762a1bSJed Brown   moments[0] = PetscRealPart(mom);
3285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetObjective(prob, 0, &f0_x));
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexComputeIntegralFEM(dm, u, &mom, user));
330c4762a1bSJed Brown   moments[1] = PetscRealPart(mom);
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetObjective(prob, 0, &f0_r2));
3325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexComputeIntegralFEM(dm, u, &mom, user));
333c4762a1bSJed Brown   moments[2] = PetscRealPart(mom);
334c4762a1bSJed Brown   PetscFunctionReturn(0);
335c4762a1bSJed Brown }
336c4762a1bSJed Brown 
337c4762a1bSJed Brown static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user)
338c4762a1bSJed Brown {
339c4762a1bSJed Brown   MPI_Comm       comm;
340c4762a1bSJed Brown   KSP            ksp;
341c4762a1bSJed Brown   Mat            M;            /* FEM mass matrix */
342c4762a1bSJed Brown   Mat            M_p;          /* Particle mass matrix */
343c4762a1bSJed Brown   Vec            f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */
344c4762a1bSJed Brown   PetscReal      pmoments[3];  /* \int f, \int x f, \int r^2 f */
345c4762a1bSJed Brown   PetscReal      fmoments[3];  /* \int \hat f, \int x \hat f, \int r^2 \hat f */
346c4762a1bSJed Brown   PetscInt       m;
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   PetscFunctionBeginUser;
3495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
3505f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(comm, &ksp));
3515f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOptionsPrefix(ksp, "ptof_"));
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &fhat));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &rhs));
354c4762a1bSJed Brown 
3555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMassMatrix(sw, dm, &M_p));
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M_p, NULL, "-M_p_view"));
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   /* make particle weight vector */
3595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* create matrix RHS vector */
3625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(M_p, f, rhs));
3635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
3645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) rhs,"rhs"));
3655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(rhs, NULL, "-rhs_view"));
366c4762a1bSJed Brown 
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(dm, &M));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
3695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M, NULL, "-M_view"));
3705f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(ksp, M, M));
3715f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(ksp));
3725f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(ksp, rhs, fhat));
3735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fhat,"fhat"));
3745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(fhat, NULL, "-fhat_view"));
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   /* Check moments of field */
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(computeParticleMoments(sw, pmoments, user));
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(computeFEMMoments(dm, fhat, fmoments, user));
3795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
380c4762a1bSJed Brown   for (m = 0; m < 3; ++m) {
3812c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol,comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol);
382c4762a1bSJed Brown   }
383c4762a1bSJed Brown 
3845f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&ksp));
3855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M));
3865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M_p));
3875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &fhat));
3885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &rhs));
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   PetscFunctionReturn(0);
391c4762a1bSJed Brown }
392c4762a1bSJed Brown 
393c4762a1bSJed Brown static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user)
394c4762a1bSJed Brown {
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   MPI_Comm       comm;
397c4762a1bSJed Brown   KSP            ksp;
398c4762a1bSJed Brown   Mat            M;            /* FEM mass matrix */
399832e2b15SMatthew G. Knepley   Mat            M_p, PM_p;    /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */
400c4762a1bSJed Brown   Vec            f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */
401c4762a1bSJed Brown   PetscReal      pmoments[3];  /* \int f, \int x f, \int r^2 f */
402c4762a1bSJed Brown   PetscReal      fmoments[3];  /* \int \hat f, \int x \hat f, \int r^2 \hat f */
403c4762a1bSJed Brown   PetscInt       m;
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   PetscFunctionBeginUser;
4065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
407c4762a1bSJed Brown 
4085f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(comm, &ksp));
4095f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOptionsPrefix(ksp, "ftop_"));
4105f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(ksp));
411c4762a1bSJed Brown 
4125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &fhat));
4135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &rhs));
414c4762a1bSJed Brown 
4155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMassMatrix(sw, dm, &M_p));
4165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M_p, NULL, "-M_p_view"));
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   /* make particle weight vector */
4195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
420c4762a1bSJed Brown 
421c4762a1bSJed Brown   /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */
4225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) rhs,"rhs"));
4235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(rhs, NULL, "-rhs_view"));
4245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(dm, &M));
4255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
4265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M, NULL, "-M_view"));
4275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(M, fhat, rhs));
4285f80ce2aSJacob Faibussowitsch   if (user->useBlockDiagPrec) CHKERRQ(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
4295f80ce2aSJacob Faibussowitsch   else                        {CHKERRQ(PetscObjectReference((PetscObject) M_p)); PM_p = M_p;}
430c4762a1bSJed Brown 
4315f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(ksp, M_p, PM_p));
4325f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolveTranspose(ksp, rhs, f));
4335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fhat,"fhat"));
4345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(fhat, NULL, "-fhat_view"));
435c4762a1bSJed Brown 
4365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   /* Check moments */
4395f80ce2aSJacob Faibussowitsch   CHKERRQ(computeParticleMoments(sw, pmoments, user));
4405f80ce2aSJacob Faibussowitsch   CHKERRQ(computeFEMMoments(dm, fhat, fmoments, user));
4415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
442c4762a1bSJed Brown   for (m = 0; m < 3; ++m) {
4432c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol,comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol);
444c4762a1bSJed Brown   }
4455f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&ksp));
4465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M));
4475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M_p));
4485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&PM_p));
4495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &fhat));
4505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &rhs));
451c4762a1bSJed Brown   PetscFunctionReturn(0);
452c4762a1bSJed Brown }
453c4762a1bSJed Brown 
454c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */
45570a7d78aSStefano Zampini static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC)
45670a7d78aSStefano Zampini {
457c4762a1bSJed Brown   DM_Plex         *mesh  = (DM_Plex *) dm->data;
458c4762a1bSJed Brown   PetscInt         debug = mesh->printFEM;
459c4762a1bSJed Brown   DM               dmC;
460c4762a1bSJed Brown   PetscSection     section;
461c4762a1bSJed Brown   PetscQuadrature  quad = NULL;
462c4762a1bSJed Brown   PetscScalar     *interpolant, *gradsum;
463c4762a1bSJed Brown   PetscFEGeom      fegeom;
464c4762a1bSJed Brown   PetscReal       *coords;
465c4762a1bSJed Brown   const PetscReal *quadPoints, *quadWeights;
466c4762a1bSJed Brown   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   PetscFunctionBegin;
4695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetDM(locC, &dmC));
4705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(locC, 0.0));
4715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
4725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dm, &coordDim));
473c4762a1bSJed Brown   fegeom.dimEmbed = coordDim;
4745f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm, &section));
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetNumFields(section, &numFields));
476c4762a1bSJed Brown   for (field = 0; field < numFields; ++field) {
477c4762a1bSJed Brown     PetscObject  obj;
478c4762a1bSJed Brown     PetscClassId id;
479c4762a1bSJed Brown     PetscInt     Nc;
480c4762a1bSJed Brown 
4815f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetField(dm, field, NULL, &obj));
4825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetClassId(obj, &id));
483c4762a1bSJed Brown     if (id == PETSCFE_CLASSID) {
484c4762a1bSJed Brown       PetscFE fe = (PetscFE) obj;
485c4762a1bSJed Brown 
4865f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFEGetQuadrature(fe, &quad));
4875f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFEGetNumComponents(fe, &Nc));
488c4762a1bSJed Brown     } else if (id == PETSCFV_CLASSID) {
489c4762a1bSJed Brown       PetscFV fv = (PetscFV) obj;
490c4762a1bSJed Brown 
4915f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFVGetQuadrature(fv, &quad));
4925f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFVGetNumComponents(fv, &Nc));
49398921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
494c4762a1bSJed Brown     numComponents += Nc;
495c4762a1bSJed Brown   }
4965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
4972c71b3e2SJacob Faibussowitsch   PetscCheckFalse((qNc != 1) && (qNc != numComponents),PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents);
4985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc6(coordDim*numComponents*2,&gradsum,coordDim*numComponents,&interpolant,coordDim*Nq,&coords,Nq,&fegeom.detJ,coordDim*coordDim*Nq,&fegeom.J,coordDim*coordDim*Nq,&fegeom.invJ));
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
5005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
501c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
502c4762a1bSJed Brown     PetscInt *star = NULL;
503c4762a1bSJed Brown     PetscInt starSize, st, d, fc;
504c4762a1bSJed Brown 
5055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArrayzero(gradsum, coordDim*numComponents));
5065f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
507c4762a1bSJed Brown     for (st = 0; st < starSize*2; st += 2) {
508c4762a1bSJed Brown       const PetscInt cell = star[st];
509c4762a1bSJed Brown       PetscScalar    *grad = &gradsum[coordDim*numComponents];
510c4762a1bSJed Brown       PetscScalar    *x    = NULL;
511c4762a1bSJed Brown 
512c4762a1bSJed Brown       if ((cell < cStart) || (cell >= cEnd)) continue;
5135f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
5145f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
515c4762a1bSJed Brown       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
516c4762a1bSJed Brown         PetscObject  obj;
517c4762a1bSJed Brown         PetscClassId id;
518c4762a1bSJed Brown         PetscInt     Nb, Nc, q, qc = 0;
519c4762a1bSJed Brown 
5205f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArrayzero(grad, coordDim*numComponents));
5215f80ce2aSJacob Faibussowitsch         CHKERRQ(DMGetField(dm, field, NULL, &obj));
5225f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscObjectGetClassId(obj, &id));
5235f80ce2aSJacob Faibussowitsch         if (id == PETSCFE_CLASSID)      {CHKERRQ(PetscFEGetNumComponents((PetscFE) obj, &Nc));CHKERRQ(PetscFEGetDimension((PetscFE) obj, &Nb));}
5245f80ce2aSJacob Faibussowitsch         else if (id == PETSCFV_CLASSID) {CHKERRQ(PetscFVGetNumComponents((PetscFV) obj, &Nc));Nb = 1;}
52598921bdaSJacob Faibussowitsch         else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
526c4762a1bSJed Brown         for (q = 0; q < Nq; ++q) {
5275f80ce2aSJacob Faibussowitsch           PetscCheck(fegeom.detJ[q] > 0.0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", (double)fegeom.detJ[q], cell, q);
5285f80ce2aSJacob Faibussowitsch           if (id == PETSCFE_CLASSID)      CHKERRQ(PetscFEInterpolateGradient_Static((PetscFE) obj, 1, &x[fieldOffset], &fegeom, q, interpolant));
52998921bdaSJacob Faibussowitsch           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
530c4762a1bSJed Brown           for (fc = 0; fc < Nc; ++fc) {
531c4762a1bSJed Brown             const PetscReal wt = quadWeights[q*qNc+qc+fc];
532c4762a1bSJed Brown 
533c4762a1bSJed Brown             for (d = 0; d < coordDim; ++d) grad[fc*coordDim+d] += interpolant[fc*dim+d]*wt*fegeom.detJ[q];
534c4762a1bSJed Brown           }
535c4762a1bSJed Brown         }
536c4762a1bSJed Brown         fieldOffset += Nb;
537c4762a1bSJed Brown       }
5385f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
539c4762a1bSJed Brown       for (fc = 0; fc < numComponents; ++fc) {
540c4762a1bSJed Brown         for (d = 0; d < coordDim; ++d) {
541c4762a1bSJed Brown           gradsum[fc*coordDim+d] += grad[fc*coordDim+d];
542c4762a1bSJed Brown         }
543c4762a1bSJed Brown       }
544c4762a1bSJed Brown       if (debug) {
5455f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Cell %D gradient: [", cell));
546c4762a1bSJed Brown         for (fc = 0; fc < numComponents; ++fc) {
547c4762a1bSJed Brown           for (d = 0; d < coordDim; ++d) {
5485f80ce2aSJacob Faibussowitsch             if (fc || d > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, ", "));
5495f80ce2aSJacob Faibussowitsch             CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc*coordDim+d])));
550c4762a1bSJed Brown           }
551c4762a1bSJed Brown         }
5525f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "]\n"));
553c4762a1bSJed Brown       }
554c4762a1bSJed Brown     }
5555f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
5565f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
557c4762a1bSJed Brown   }
5585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree6(gradsum,interpolant,coords,fegeom.detJ,fegeom.J,fegeom.invJ));
559c4762a1bSJed Brown   PetscFunctionReturn(0);
560c4762a1bSJed Brown }
561c4762a1bSJed Brown 
562c4762a1bSJed Brown static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user)
563c4762a1bSJed Brown {
564c4762a1bSJed Brown 
565c4762a1bSJed Brown   MPI_Comm       comm;
566c4762a1bSJed Brown   KSP            ksp;
567c4762a1bSJed Brown   Mat            M;                   /* FEM mass matrix */
568c4762a1bSJed Brown   Mat            M_p;                 /* Particle mass matrix */
569c4762a1bSJed Brown   Vec            f, rhs, fhat, grad;  /* Particle field f, \int phi_i f, FEM field */
570c4762a1bSJed Brown   PetscReal      pmoments[3];         /* \int f, \int x f, \int r^2 f */
571c4762a1bSJed Brown   PetscReal      fmoments[3];         /* \int \hat f, \int x \hat f, \int r^2 \hat f */
572c4762a1bSJed Brown   PetscInt       m;
573c4762a1bSJed Brown 
574c4762a1bSJed Brown   PetscFunctionBeginUser;
5755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
5765f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(comm, &ksp));
5775f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOptionsPrefix(ksp, "ptof_"));
5785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &fhat));
5795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &rhs));
5805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &grad));
581c4762a1bSJed Brown 
5825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMassMatrix(sw, dm, &M_p));
5835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M_p, NULL, "-M_p_view"));
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   /* make particle weight vector */
5865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   /* create matrix RHS vector */
5895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(M_p, f, rhs));
5905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
5915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) rhs,"rhs"));
5925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(rhs, NULL, "-rhs_view"));
593c4762a1bSJed Brown 
5945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(dm, &M));
5955f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
596c4762a1bSJed Brown 
5975f80ce2aSJacob Faibussowitsch   CHKERRQ(InterpolateGradient(dm, fhat, grad));
598c4762a1bSJed Brown 
5995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(M, NULL, "-M_view"));
6005f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(ksp, M, M));
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(ksp));
6025f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(ksp, rhs, grad));
6035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fhat,"fhat"));
6045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(fhat, NULL, "-fhat_view"));
605c4762a1bSJed Brown 
606c4762a1bSJed Brown   /* Check moments of field */
6075f80ce2aSJacob Faibussowitsch   CHKERRQ(computeParticleMoments(sw, pmoments, user));
6085f80ce2aSJacob Faibussowitsch   CHKERRQ(computeFEMMoments(dm, grad, fmoments, user));
6095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
610c4762a1bSJed Brown   for (m = 0; m < 3; ++m) {
6112c71b3e2SJacob Faibussowitsch     PetscCheckFalse(PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]) > user->momentTol,comm, PETSC_ERR_ARG_WRONG, "Moment %D error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m])/fmoments[m]), user->momentTol);
612c4762a1bSJed Brown   }
613c4762a1bSJed Brown 
6145f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&ksp));
6155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M));
6165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&M_p));
6175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &fhat));
6185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &rhs));
6195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &grad));
620c4762a1bSJed Brown 
621c4762a1bSJed Brown   PetscFunctionReturn(0);
622c4762a1bSJed Brown }
623c4762a1bSJed Brown 
6245627991aSBarry Smith int main (int argc, char *argv[])
6255627991aSBarry Smith {
626c4762a1bSJed Brown   MPI_Comm       comm;
627c4762a1bSJed Brown   DM             dm, sw;
628c4762a1bSJed Brown   AppCtx         user;
629c4762a1bSJed Brown 
630*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
631c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
6325f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(comm, &user));
6335f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &dm, &user));
6345f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateFEM(dm, &user));
6355f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateParticles(dm, &sw, &user));
6365f80ce2aSJacob Faibussowitsch   CHKERRQ(TestL2ProjectionParticlesToField(dm, sw, &user));
6375f80ce2aSJacob Faibussowitsch   CHKERRQ(TestL2ProjectionFieldToParticles(dm, sw, &user));
6385f80ce2aSJacob Faibussowitsch   CHKERRQ(TestFieldGradientProjection(dm, sw, &user));
6395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
6405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&sw));
641*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
642*b122ec5aSJacob Faibussowitsch   return 0;
643c4762a1bSJed Brown }
644c4762a1bSJed Brown 
645c4762a1bSJed Brown /*TEST
646c4762a1bSJed Brown 
647832e2b15SMatthew G. Knepley   # Swarm does not handle complex or quad
648832e2b15SMatthew G. Knepley   build:
649832e2b15SMatthew G. Knepley     requires: !complex double
650832e2b15SMatthew G. Knepley 
651c4762a1bSJed Brown   test:
652c4762a1bSJed Brown     suffix: proj_tri_0
653832e2b15SMatthew G. Knepley     requires: triangle
654832e2b15SMatthew 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
655c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   test:
658c4762a1bSJed Brown     suffix: proj_tri_2_faces
659832e2b15SMatthew G. Knepley     requires: triangle
660832e2b15SMatthew 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
661c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
662c4762a1bSJed Brown 
663c4762a1bSJed Brown   test:
664c4762a1bSJed Brown     suffix: proj_quad_0
665832e2b15SMatthew G. Knepley     requires: triangle
66630602db0SMatthew 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
667c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
668c4762a1bSJed Brown 
669c4762a1bSJed Brown   test:
670c4762a1bSJed Brown     suffix: proj_quad_2_faces
671832e2b15SMatthew G. Knepley     requires: triangle
67230602db0SMatthew 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
673c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
674c4762a1bSJed Brown 
675c4762a1bSJed Brown   test:
676c4762a1bSJed Brown     suffix: proj_tri_5P
677832e2b15SMatthew G. Knepley     requires: triangle
678832e2b15SMatthew 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
679c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
680c4762a1bSJed Brown 
681c4762a1bSJed Brown   test:
682c4762a1bSJed Brown     suffix: proj_quad_5P
683832e2b15SMatthew G. Knepley     requires: triangle
68430602db0SMatthew 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
685c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
686c4762a1bSJed Brown 
687c4762a1bSJed Brown   test:
688c4762a1bSJed Brown     suffix: proj_tri_mdx
689832e2b15SMatthew G. Knepley     requires: triangle
690832e2b15SMatthew 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
691c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
692c4762a1bSJed Brown 
693c4762a1bSJed Brown   test:
694c4762a1bSJed Brown     suffix: proj_tri_mdx_5P
695832e2b15SMatthew G. Knepley     requires: triangle
696832e2b15SMatthew 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
697c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
698c4762a1bSJed Brown 
699c4762a1bSJed Brown   test:
700c4762a1bSJed Brown     suffix: proj_tri_3d
701832e2b15SMatthew G. Knepley     requires: ctetgen
70230602db0SMatthew 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
703c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   test:
706c4762a1bSJed Brown     suffix: proj_tri_3d_2_faces
707832e2b15SMatthew G. Knepley     requires: ctetgen
70830602db0SMatthew 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
709c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
710c4762a1bSJed Brown 
711c4762a1bSJed Brown   test:
712c4762a1bSJed Brown     suffix: proj_tri_3d_5P
713832e2b15SMatthew G. Knepley     requires: ctetgen
71430602db0SMatthew 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
715c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
716c4762a1bSJed Brown 
717c4762a1bSJed Brown   test:
718c4762a1bSJed Brown     suffix: proj_tri_3d_mdx
719832e2b15SMatthew G. Knepley     requires: ctetgen
72030602db0SMatthew 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
721c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
722c4762a1bSJed Brown 
723c4762a1bSJed Brown   test:
724c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_5P
725832e2b15SMatthew G. Knepley     requires: ctetgen
72630602db0SMatthew 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
727c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
728c4762a1bSJed Brown 
729c4762a1bSJed Brown   test:
730c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_2_faces
731832e2b15SMatthew G. Knepley     requires: ctetgen
73230602db0SMatthew 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
733c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
734c4762a1bSJed Brown 
735c4762a1bSJed Brown   test:
736c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_5P_2_faces
737832e2b15SMatthew G. Knepley     requires: ctetgen
73830602db0SMatthew 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
739c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
740c4762a1bSJed Brown 
741832e2b15SMatthew G. Knepley   test:
742832e2b15SMatthew G. Knepley     suffix: proj_quad_lsqr_scale
743832e2b15SMatthew G. Knepley     requires: !complex
74430602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
745832e2b15SMatthew G. Knepley       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
746832e2b15SMatthew G. Knepley       -particlesPerCell 200 \
747832e2b15SMatthew G. Knepley       -ptof_pc_type lu  \
748832e2b15SMatthew G. Knepley       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none
749832e2b15SMatthew G. Knepley 
750832e2b15SMatthew G. Knepley   test:
751832e2b15SMatthew G. Knepley     suffix: proj_quad_lsqr_prec_scale
752832e2b15SMatthew G. Knepley     requires: !complex
75330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
754832e2b15SMatthew G. Knepley       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
755832e2b15SMatthew G. Knepley       -particlesPerCell 200 \
756832e2b15SMatthew G. Knepley       -ptof_pc_type lu  \
757832e2b15SMatthew G. Knepley       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec
758832e2b15SMatthew G. Knepley 
759c4762a1bSJed Brown TEST*/
760