xref: /petsc/src/dm/impls/swarm/tests/ex2.c (revision 09b68a49ed2854d1e4985cc2aa6af33c7c4e69b3)
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 
linear(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,void * a_ctx)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 
x2_x4(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,void * a_ctx)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 
sinx(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,void * a_ctx)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 
ProcessOptions(MPI_Comm comm,AppCtx * options)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();
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
PerturbVertices(DM dm,AppCtx * user)88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PerturbVertices(DM dm, AppCtx *user)
89d71ae5a4SJacob Faibussowitsch {
90c4762a1bSJed Brown   PetscRandom  rnd;
91c4762a1bSJed Brown   PetscReal    interval = user->meshRelDx;
92c4762a1bSJed Brown   Vec          coordinates;
93c4762a1bSJed Brown   PetscScalar *coords;
94832e2b15SMatthew G. Knepley   PetscReal   *hh, low[3], high[3];
95832e2b15SMatthew G. Knepley   PetscInt     d, cdim, cEnd, N, p, bs;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   PetscFunctionBeginUser;
989566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dm, low, high));
999566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd));
1009566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
1019566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -interval, interval));
1029566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
1039566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1049566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
1059566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(cdim, &hh));
106832e2b15SMatthew G. Knepley   for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim);
1079566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(coordinates, &N));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(coordinates, &bs));
10963a3b9bcSJacob Faibussowitsch   PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim);
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
111c4762a1bSJed Brown   for (p = 0; p < N; p += cdim) {
112c4762a1bSJed Brown     PetscScalar *coord = &coords[p], value;
113c4762a1bSJed Brown 
114c4762a1bSJed Brown     for (d = 0; d < cdim; ++d) {
1159566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rnd, &value));
116832e2b15SMatthew G. Knepley       coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d])));
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown   }
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
1209566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
1219566063dSJacob Faibussowitsch   PetscCall(PetscFree(hh));
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * user)125d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
126d71ae5a4SJacob Faibussowitsch {
127832e2b15SMatthew G. Knepley   PetscReal low[3], high[3];
128832e2b15SMatthew G. Knepley   PetscInt  cdim, d;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   PetscFunctionBeginUser;
1319566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1329566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1339566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1349566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
13530602db0SMatthew G. Knepley 
1369566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(*dm, &cdim));
1379566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(*dm, low, high));
138832e2b15SMatthew G. Knepley   for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d];
1399566063dSJacob Faibussowitsch   PetscCall(PerturbVertices(*dm, user));
1403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
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[])143d71ae5a4SJacob 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[])
144d71ae5a4SJacob Faibussowitsch {
145c4762a1bSJed Brown   g0[0] = 1.0;
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
CreateFEM(DM dm,AppCtx * user)148d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
149d71ae5a4SJacob Faibussowitsch {
150c4762a1bSJed Brown   PetscFE        fe;
151c4762a1bSJed Brown   PetscDS        ds;
152832e2b15SMatthew G. Knepley   DMPolytopeType ct;
153832e2b15SMatthew G. Knepley   PetscInt       dim, cStart;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   PetscFunctionBeginUser;
1569566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1579566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
159fc3eb83cSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PetscObjectComm((PetscObject)dm), dim, 1, ct, NULL, -1, &fe));
1609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
1619566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1629566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
164c4762a1bSJed Brown   /* Setup to form mass matrix */
1659566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1669566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL));
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
CreateParticles(DM dm,DM * sw,AppCtx * user)170d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
171d71ae5a4SJacob Faibussowitsch {
17219307e5cSMatthew G. Knepley   DMSwarmCellDM  celldm;
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;
17919307e5cSMatthew G. Knepley   PetscInt      *swarm_cellid;
18019307e5cSMatthew G. Knepley   PetscInt       Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d, Nfc;
18119307e5cSMatthew G. Knepley   const char   **coordFields, *cellid;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   PetscFunctionBeginUser;
1849566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1869566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
187832e2b15SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
1889566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1899566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1909566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
191c4762a1bSJed Brown 
1929566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd));
1939566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0));
1949566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
1959566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rndp));
1969566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rndp, -interval, interval));
1979566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rndp));
198c4762a1bSJed Brown 
1999566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
2009566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
2019566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
2029566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
2039566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell));
2049566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0));
2059566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
2069072cb8bSMatthew G. Knepley 
2079072cb8bSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(*sw, &celldm));
2089072cb8bSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
2099072cb8bSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
2109072cb8bSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
2119072cb8bSMatthew G. Knepley 
21219307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, coordFields[0], NULL, NULL, (void **)&coords));
21319307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
2149566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals));
215c4762a1bSJed Brown 
2169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
217c4762a1bSJed Brown   for (c = 0; c < Ncell; ++c) {
218c4762a1bSJed Brown     if (Np == 1) {
2199566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
22019307e5cSMatthew G. Knepley       swarm_cellid[c] = c;
221c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d];
222c4762a1bSJed Brown     } else {
2239566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
224c4762a1bSJed Brown       for (d = 0; d < dim; ++d) xi0[d] = -1.0;
225c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
226c4762a1bSJed Brown         const PetscInt n   = c * Np + p;
227c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
228c4762a1bSJed Brown 
22919307e5cSMatthew G. Knepley         swarm_cellid[n] = c;
2309371c9d4SSatish Balay         for (d = 0; d < dim; ++d) {
2319371c9d4SSatish Balay           PetscCall(PetscRandomGetValue(rnd, &value));
2329371c9d4SSatish Balay           refcoords[d] = PetscRealPart(value);
2339371c9d4SSatish Balay           sum += refcoords[d];
2349371c9d4SSatish Balay         }
2359371c9d4SSatish Balay         if (simplex && sum > 0.0)
2369371c9d4SSatish Balay           for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum;
237c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]);
238c4762a1bSJed Brown       }
239c4762a1bSJed Brown     }
240c4762a1bSJed Brown   }
2419566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
242c4762a1bSJed Brown   for (c = 0; c < Ncell; ++c) {
243c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
244c4762a1bSJed Brown       const PetscInt n = c * Np + p;
245c4762a1bSJed Brown 
2469371c9d4SSatish Balay       for (d = 0; d < dim; ++d) {
2479371c9d4SSatish Balay         PetscCall(PetscRandomGetValue(rndp, &value));
2489371c9d4SSatish Balay         coords[n * dim + d] += PetscRealPart(value);
2499371c9d4SSatish Balay       }
2503ba16761SJacob Faibussowitsch       PetscCall(user->func(dim, 0.0, &coords[n * dim], 1, &vals[c], user));
251c4762a1bSJed Brown     }
252c4762a1bSJed Brown   }
25319307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, coordFields[0], NULL, NULL, (void **)&coords));
25419307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
2559566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals));
2569566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
2579566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rndp));
2589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
259d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmVectorDefineField(*sw, "w_q"));
2609566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
2613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
CreateParticles_Shape(DM dm,DM * sw,AppCtx * user)26480aecfd0Sjosephpu static PetscErrorCode CreateParticles_Shape(DM dm, DM *sw, AppCtx *user)
26580aecfd0Sjosephpu {
26619307e5cSMatthew G. Knepley   DMSwarmCellDM    celldm;
26780aecfd0Sjosephpu   PetscDS          prob;
26880aecfd0Sjosephpu   PetscFE          fe;
26980aecfd0Sjosephpu   PetscQuadrature  quad;
27080aecfd0Sjosephpu   PetscScalar     *vals;
27180aecfd0Sjosephpu   PetscReal       *v0, *J, *invJ, detJ, *coords, *xi0;
2729072cb8bSMatthew G. Knepley   PetscInt        *swarm_cellid;
27380aecfd0Sjosephpu   const PetscReal *qpoints, *qweights;
27419307e5cSMatthew G. Knepley   PetscInt         cStart, cEnd, c, Nq, q, dim, Nfc;
27519307e5cSMatthew G. Knepley   const char     **coordFields, *cellid;
27680aecfd0Sjosephpu 
27780aecfd0Sjosephpu   PetscFunctionBeginUser;
27880aecfd0Sjosephpu   PetscCall(DMGetDimension(dm, &dim));
27980aecfd0Sjosephpu   PetscCall(DMGetDS(dm, &prob));
28080aecfd0Sjosephpu   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
28180aecfd0Sjosephpu   PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe));
28280aecfd0Sjosephpu   PetscCall(PetscFEGetQuadrature(fe, &quad));
28380aecfd0Sjosephpu   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, &qweights));
28480aecfd0Sjosephpu   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
28580aecfd0Sjosephpu   PetscCall(DMSetType(*sw, DMSWARM));
28680aecfd0Sjosephpu   PetscCall(DMSetDimension(*sw, dim));
28780aecfd0Sjosephpu   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
28880aecfd0Sjosephpu   PetscCall(DMSwarmSetCellDM(*sw, dm));
2899072cb8bSMatthew G. Knepley   PetscCall(DMSwarmGetCellDMActive(*sw, &celldm));
29019307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCellID(celldm, &cellid));
29119307e5cSMatthew G. Knepley   PetscCall(DMSwarmCellDMGetCoordinateFields(celldm, &Nfc, &coordFields));
29219307e5cSMatthew G. Knepley   PetscCheck(Nfc == 1, PetscObjectComm((PetscObject)sw), PETSC_ERR_SUP, "We only support a single coordinate field right now, not %" PetscInt_FMT, Nfc);
29380aecfd0Sjosephpu   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_REAL));
29480aecfd0Sjosephpu   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
29580aecfd0Sjosephpu   PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Nq, 0));
29680aecfd0Sjosephpu   PetscCall(DMSetFromOptions(*sw));
29780aecfd0Sjosephpu   PetscCall(PetscMalloc4(dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ));
29880aecfd0Sjosephpu   for (c = 0; c < dim; c++) xi0[c] = -1.;
29919307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, coordFields[0], NULL, NULL, (void **)&coords));
30019307e5cSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
30180aecfd0Sjosephpu   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals));
30280aecfd0Sjosephpu   for (c = cStart; c < cEnd; ++c) {
30380aecfd0Sjosephpu     for (q = 0; q < Nq; ++q) {
30480aecfd0Sjosephpu       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ));
30519307e5cSMatthew G. Knepley       swarm_cellid[c * Nq + q] = c;
30680aecfd0Sjosephpu       CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], &coords[(c * Nq + q) * dim]);
30780aecfd0Sjosephpu       PetscCall(linear(dim, 0.0, &coords[(c * Nq + q) * dim], 1, &vals[c * Nq + q], (void *)user));
30880aecfd0Sjosephpu       vals[c * Nq + q] *= qweights[q] * detJ;
30980aecfd0Sjosephpu     }
31080aecfd0Sjosephpu   }
31119307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, coordFields[0], NULL, NULL, (void **)&coords));
31219307e5cSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, cellid, NULL, NULL, (void **)&swarm_cellid));
31380aecfd0Sjosephpu   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals));
31480aecfd0Sjosephpu   PetscCall(PetscFree4(xi0, v0, J, invJ));
31580aecfd0Sjosephpu   PetscCall(DMSwarmMigrate(*sw, PETSC_FALSE));
31680aecfd0Sjosephpu   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
317d52c2f21SMatthew G. Knepley   PetscCall(DMSwarmVectorDefineField(*sw, "w_q"));
31880aecfd0Sjosephpu   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
31980aecfd0Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
32080aecfd0Sjosephpu }
32180aecfd0Sjosephpu 
computeParticleMoments(DM sw,PetscReal moments[3],AppCtx * user)322d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user)
323d71ae5a4SJacob Faibussowitsch {
324c4762a1bSJed Brown   DM                 dm;
325c4762a1bSJed Brown   const PetscReal   *coords;
326c4762a1bSJed Brown   const PetscScalar *w;
327c4762a1bSJed Brown   PetscReal          mom[3] = {0.0, 0.0, 0.0};
328c4762a1bSJed Brown   PetscInt           cell, cStart, cEnd, dim;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBeginUser;
3319566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sw, &dim));
3329566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetCellDM(sw, &dm));
3339566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3349566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(sw));
3359566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
3369566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w));
337c4762a1bSJed Brown   for (cell = cStart; cell < cEnd; ++cell) {
338c4762a1bSJed Brown     PetscInt *pidx;
339c4762a1bSJed Brown     PetscInt  Np, p, d;
340c4762a1bSJed Brown 
3419566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
342c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
343c4762a1bSJed Brown       const PetscInt   idx = pidx[p];
344c4762a1bSJed Brown       const PetscReal *c   = &coords[idx * dim];
345c4762a1bSJed Brown 
346c4762a1bSJed Brown       mom[0] += PetscRealPart(w[idx]);
347c4762a1bSJed Brown       mom[1] += PetscRealPart(w[idx]) * c[0];
348c4762a1bSJed Brown       for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d];
349c4762a1bSJed Brown     }
350fe11354eSMatthew G. Knepley     PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Np, &pidx));
351c4762a1bSJed Brown   }
3529566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
3539566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w));
3549566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(sw));
355462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
3563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
357c4762a1bSJed Brown }
358c4762a1bSJed Brown 
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[])359d71ae5a4SJacob 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[])
360d71ae5a4SJacob Faibussowitsch {
361c4762a1bSJed Brown   f0[0] = u[0];
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
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[])364d71ae5a4SJacob 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[])
365d71ae5a4SJacob Faibussowitsch {
366c4762a1bSJed Brown   f0[0] = x[0] * u[0];
367c4762a1bSJed Brown }
368c4762a1bSJed Brown 
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[])369d71ae5a4SJacob 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[])
370d71ae5a4SJacob Faibussowitsch {
371c4762a1bSJed Brown   PetscInt d;
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   f0[0] = 0.0;
374c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0];
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
computeFEMMoments(DM dm,Vec u,PetscReal moments[3],AppCtx * user)377d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
378d71ae5a4SJacob Faibussowitsch {
379c4762a1bSJed Brown   PetscDS     prob;
380c4762a1bSJed Brown   PetscScalar mom;
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   PetscFunctionBeginUser;
3839566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
3849566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_1));
3859566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
386c4762a1bSJed Brown   moments[0] = PetscRealPart(mom);
3879566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_x));
3889566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
389c4762a1bSJed Brown   moments[1] = PetscRealPart(mom);
3909566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_r2));
3919566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
392c4762a1bSJed Brown   moments[2] = PetscRealPart(mom);
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
394c4762a1bSJed Brown }
395c4762a1bSJed Brown 
TestL2ProjectionParticlesToField(DM dm,DM sw,Vec fhat,AppCtx * user)396fc3eb83cSMatthew G. Knepley static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, Vec fhat, AppCtx *user)
397d71ae5a4SJacob Faibussowitsch {
398fc3eb83cSMatthew G. Knepley   const char *fieldnames[1] = {"w_q"};
399fc3eb83cSMatthew G. Knepley   Vec         fields[1]     = {fhat};
400fc3eb83cSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
401fc3eb83cSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   PetscFunctionBeginUser;
404953494dbSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, dm, 1, fieldnames, fields, SCATTER_FORWARD));
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   /* Check moments of field */
4079566063dSJacob Faibussowitsch   PetscCall(computeParticleMoments(sw, pmoments, user));
4089566063dSJacob Faibussowitsch   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
409fc3eb83cSMatthew 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]));
410fc3eb83cSMatthew G. Knepley   for (PetscInt m = 0; m < 3; ++m) {
411fc3eb83cSMatthew 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]),
412fc3eb83cSMatthew G. Knepley                user->momentTol);
413c4762a1bSJed Brown   }
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415c4762a1bSJed Brown }
416c4762a1bSJed Brown 
TestL2ProjectionFieldToParticles(DM dm,DM sw,Vec fhat,AppCtx * user)417fc3eb83cSMatthew G. Knepley static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, Vec fhat, AppCtx *user)
418d71ae5a4SJacob Faibussowitsch {
419fc3eb83cSMatthew G. Knepley   const char *fieldnames[1] = {"w_q"};
420fc3eb83cSMatthew G. Knepley   Vec         fields[1]     = {fhat};
421fc3eb83cSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
422fc3eb83cSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   PetscFunctionBeginUser;
425953494dbSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, dm, 1, fieldnames, fields, SCATTER_REVERSE));
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /* Check moments */
4289566063dSJacob Faibussowitsch   PetscCall(computeParticleMoments(sw, pmoments, user));
4299566063dSJacob Faibussowitsch   PetscCall(computeFEMMoments(dm, fhat, fmoments, user));
430fc3eb83cSMatthew 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]));
431fc3eb83cSMatthew G. Knepley   for (PetscInt m = 0; m < 3; ++m) {
432fc3eb83cSMatthew 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]),
433fc3eb83cSMatthew G. Knepley                user->momentTol);
434c4762a1bSJed Brown   }
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
436c4762a1bSJed Brown }
437c4762a1bSJed Brown 
438c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */
InterpolateGradient(DM dm,Vec locX,Vec locC)439d71ae5a4SJacob Faibussowitsch static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC)
440d71ae5a4SJacob Faibussowitsch {
441c4762a1bSJed Brown   DM_Plex         *mesh  = (DM_Plex *)dm->data;
442c4762a1bSJed Brown   PetscInt         debug = mesh->printFEM;
443c4762a1bSJed Brown   DM               dmC;
444c4762a1bSJed Brown   PetscSection     section;
445c4762a1bSJed Brown   PetscQuadrature  quad = NULL;
446c4762a1bSJed Brown   PetscScalar     *interpolant, *gradsum;
447c4762a1bSJed Brown   PetscFEGeom      fegeom;
448c4762a1bSJed Brown   PetscReal       *coords;
449c4762a1bSJed Brown   const PetscReal *quadPoints, *quadWeights;
450c4762a1bSJed Brown   PetscInt         dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset;
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   PetscFunctionBegin;
4539566063dSJacob Faibussowitsch   PetscCall(VecGetDM(locC, &dmC));
4549566063dSJacob Faibussowitsch   PetscCall(VecSet(locC, 0.0));
4559566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4569566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &coordDim));
457c4762a1bSJed Brown   fegeom.dimEmbed = coordDim;
4589566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
4599566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &numFields));
460c4762a1bSJed Brown   for (field = 0; field < numFields; ++field) {
461c4762a1bSJed Brown     PetscObject  obj;
462c4762a1bSJed Brown     PetscClassId id;
463c4762a1bSJed Brown     PetscInt     Nc;
464c4762a1bSJed Brown 
4659566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, field, NULL, &obj));
4669566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetClassId(obj, &id));
467c4762a1bSJed Brown     if (id == PETSCFE_CLASSID) {
468c4762a1bSJed Brown       PetscFE fe = (PetscFE)obj;
469c4762a1bSJed Brown 
4709566063dSJacob Faibussowitsch       PetscCall(PetscFEGetQuadrature(fe, &quad));
4719566063dSJacob Faibussowitsch       PetscCall(PetscFEGetNumComponents(fe, &Nc));
472c4762a1bSJed Brown     } else if (id == PETSCFV_CLASSID) {
473c4762a1bSJed Brown       PetscFV fv = (PetscFV)obj;
474c4762a1bSJed Brown 
4759566063dSJacob Faibussowitsch       PetscCall(PetscFVGetQuadrature(fv, &quad));
4769566063dSJacob Faibussowitsch       PetscCall(PetscFVGetNumComponents(fv, &Nc));
47763a3b9bcSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
478c4762a1bSJed Brown     numComponents += Nc;
479c4762a1bSJed Brown   }
4809566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights));
48163a3b9bcSJacob Faibussowitsch   PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents);
4829566063dSJacob 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));
4839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4849566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
485c4762a1bSJed Brown   for (v = vStart; v < vEnd; ++v) {
486c4762a1bSJed Brown     PetscInt *star = NULL;
487c4762a1bSJed Brown     PetscInt  starSize, st, d, fc;
488c4762a1bSJed Brown 
4899566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(gradsum, coordDim * numComponents));
4909566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
491c4762a1bSJed Brown     for (st = 0; st < starSize * 2; st += 2) {
492c4762a1bSJed Brown       const PetscInt cell = star[st];
493c4762a1bSJed Brown       PetscScalar   *grad = &gradsum[coordDim * numComponents];
494c4762a1bSJed Brown       PetscScalar   *x    = NULL;
495c4762a1bSJed Brown 
496c4762a1bSJed Brown       if ((cell < cStart) || (cell >= cEnd)) continue;
4979566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ));
4989566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x));
499c4762a1bSJed Brown       for (field = 0, fieldOffset = 0; field < numFields; ++field) {
500c4762a1bSJed Brown         PetscObject  obj;
501c4762a1bSJed Brown         PetscClassId id;
502c4762a1bSJed Brown         PetscInt     Nb, Nc, q, qc = 0;
503c4762a1bSJed Brown 
5049566063dSJacob Faibussowitsch         PetscCall(PetscArrayzero(grad, coordDim * numComponents));
5059566063dSJacob Faibussowitsch         PetscCall(DMGetField(dm, field, NULL, &obj));
5069566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetClassId(obj, &id));
5079371c9d4SSatish Balay         if (id == PETSCFE_CLASSID) {
5089371c9d4SSatish Balay           PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc));
5099371c9d4SSatish Balay           PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb));
5109371c9d4SSatish Balay         } else if (id == PETSCFV_CLASSID) {
5119371c9d4SSatish Balay           PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc));
5129371c9d4SSatish Balay           Nb = 1;
5139371c9d4SSatish Balay         } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
514c4762a1bSJed Brown         for (q = 0; q < Nq; ++q) {
51563a3b9bcSJacob 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);
5169566063dSJacob Faibussowitsch           if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant));
51763a3b9bcSJacob Faibussowitsch           else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field);
518c4762a1bSJed Brown           for (fc = 0; fc < Nc; ++fc) {
519c4762a1bSJed Brown             const PetscReal wt = quadWeights[q * qNc + qc + fc];
520c4762a1bSJed Brown 
521c4762a1bSJed Brown             for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q];
522c4762a1bSJed Brown           }
523c4762a1bSJed Brown         }
524c4762a1bSJed Brown         fieldOffset += Nb;
525c4762a1bSJed Brown       }
5269566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x));
527c4762a1bSJed Brown       for (fc = 0; fc < numComponents; ++fc) {
528ad540459SPierre Jolivet         for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d];
529c4762a1bSJed Brown       }
530c4762a1bSJed Brown       if (debug) {
53163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell));
532c4762a1bSJed Brown         for (fc = 0; fc < numComponents; ++fc) {
533c4762a1bSJed Brown           for (d = 0; d < coordDim; ++d) {
5349566063dSJacob Faibussowitsch             if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
5359566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d])));
536c4762a1bSJed Brown           }
537c4762a1bSJed Brown         }
5389566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n"));
539c4762a1bSJed Brown       }
540c4762a1bSJed Brown     }
5419566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star));
5429566063dSJacob Faibussowitsch     PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES));
543c4762a1bSJed Brown   }
5449566063dSJacob Faibussowitsch   PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ));
5453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
546c4762a1bSJed Brown }
547c4762a1bSJed Brown 
TestFieldGradientProjection(DM dm,DM sw,AppCtx * user)548d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user)
549d71ae5a4SJacob Faibussowitsch {
550c4762a1bSJed Brown   MPI_Comm  comm;
551c4762a1bSJed Brown   KSP       ksp;
552c4762a1bSJed Brown   Mat       M;                  /* FEM mass matrix */
553c4762a1bSJed Brown   Mat       M_p;                /* Particle mass matrix */
554c4762a1bSJed Brown   Vec       f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */
555c4762a1bSJed Brown   PetscReal pmoments[3];        /* \int f, \int x f, \int r^2 f */
556c4762a1bSJed Brown   PetscReal fmoments[3];        /* \int \hat f, \int x \hat f, \int r^2 \hat f */
557c4762a1bSJed Brown   PetscInt  m;
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   PetscFunctionBeginUser;
5609566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5619566063dSJacob Faibussowitsch   PetscCall(KSPCreate(comm, &ksp));
5629566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(ksp, "ptof_"));
5639566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &fhat));
5649566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &rhs));
5659566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &grad));
566c4762a1bSJed Brown 
5679566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
5689566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view"));
569c4762a1bSJed Brown 
570c4762a1bSJed Brown   /* make particle weight vector */
5719566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   /* create matrix RHS vector */
5749566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(M_p, f, rhs));
5759566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
5769566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs"));
5779566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
578c4762a1bSJed Brown 
5799566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &M));
5809566063dSJacob Faibussowitsch   PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user));
581c4762a1bSJed Brown 
5829566063dSJacob Faibussowitsch   PetscCall(InterpolateGradient(dm, fhat, grad));
583c4762a1bSJed Brown 
5849566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M, NULL, "-M_view"));
5859566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, M, M));
5869566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
5879566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, rhs, grad));
5889566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat"));
5899566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   /* Check moments of field */
5929566063dSJacob Faibussowitsch   PetscCall(computeParticleMoments(sw, pmoments, user));
5939566063dSJacob Faibussowitsch   PetscCall(computeFEMMoments(dm, grad, fmoments, user));
5949566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
595c4762a1bSJed Brown   for (m = 0; m < 3; ++m) {
5961dca8a05SBarry 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);
597c4762a1bSJed Brown   }
598c4762a1bSJed Brown 
5999566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
6009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M));
6019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M_p));
6029566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &fhat));
6039566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &rhs));
6049566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &grad));
6053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
606c4762a1bSJed Brown }
607c4762a1bSJed Brown 
TestL2Projection_Shape(DM dm,DM sw,AppCtx * user)60880aecfd0Sjosephpu static PetscErrorCode TestL2Projection_Shape(DM dm, DM sw, AppCtx *user)
60980aecfd0Sjosephpu {
61080aecfd0Sjosephpu   PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
61180aecfd0Sjosephpu   KSP       ksp;
61280aecfd0Sjosephpu   Mat       mass;
61380aecfd0Sjosephpu   Vec       u, rhs, uproj;
61480aecfd0Sjosephpu   void     *ctxs[1];
61580aecfd0Sjosephpu   PetscReal error;
61680aecfd0Sjosephpu 
61780aecfd0Sjosephpu   PetscFunctionBeginUser;
61880aecfd0Sjosephpu   ctxs[0]  = (void *)user;
61980aecfd0Sjosephpu   funcs[0] = linear;
62080aecfd0Sjosephpu 
62180aecfd0Sjosephpu   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &u));
62280aecfd0Sjosephpu   PetscCall(VecViewFromOptions(u, NULL, "-f_view"));
62380aecfd0Sjosephpu   PetscCall(DMGetGlobalVector(dm, &rhs));
62480aecfd0Sjosephpu   PetscCall(DMCreateMassMatrix(sw, dm, &mass));
62580aecfd0Sjosephpu   PetscCall(MatMultTranspose(mass, u, rhs));
62680aecfd0Sjosephpu   PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
62780aecfd0Sjosephpu   PetscCall(MatDestroy(&mass));
62880aecfd0Sjosephpu   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &u));
62980aecfd0Sjosephpu   PetscCall(DMGetGlobalVector(dm, &uproj));
63080aecfd0Sjosephpu   PetscCall(DMCreateMatrix(dm, &mass));
63180aecfd0Sjosephpu   PetscCall(DMPlexSNESComputeJacobianFEM(dm, uproj, mass, mass, user));
63280aecfd0Sjosephpu   PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view"));
63380aecfd0Sjosephpu   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
63480aecfd0Sjosephpu   PetscCall(KSPSetOperators(ksp, mass, mass));
63580aecfd0Sjosephpu   PetscCall(KSPSetFromOptions(ksp));
63680aecfd0Sjosephpu   PetscCall(KSPSolve(ksp, rhs, uproj));
63780aecfd0Sjosephpu   PetscCall(KSPDestroy(&ksp));
63880aecfd0Sjosephpu   PetscCall(PetscObjectSetName((PetscObject)uproj, "Full Projection"));
63980aecfd0Sjosephpu   PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view"));
64080aecfd0Sjosephpu   PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, uproj, &error));
64180aecfd0Sjosephpu   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error));
64280aecfd0Sjosephpu   PetscCall(MatDestroy(&mass));
64380aecfd0Sjosephpu   PetscCall(DMRestoreGlobalVector(dm, &rhs));
64480aecfd0Sjosephpu   PetscCall(DMRestoreGlobalVector(dm, &uproj));
64580aecfd0Sjosephpu   PetscFunctionReturn(PETSC_SUCCESS);
64680aecfd0Sjosephpu }
64780aecfd0Sjosephpu 
main(int argc,char * argv[])648d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
649d71ae5a4SJacob Faibussowitsch {
650c4762a1bSJed Brown   MPI_Comm comm;
651c4762a1bSJed Brown   DM       dm, sw;
652c4762a1bSJed Brown   AppCtx   user;
653c4762a1bSJed Brown 
654327415f7SBarry Smith   PetscFunctionBeginUser;
6559566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
656c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
6579566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
6589566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
6599566063dSJacob Faibussowitsch   PetscCall(CreateFEM(dm, &user));
66080aecfd0Sjosephpu   if (!user.shape) {
661fc3eb83cSMatthew G. Knepley     Vec fhat;
662fc3eb83cSMatthew G. Knepley 
6639566063dSJacob Faibussowitsch     PetscCall(CreateParticles(dm, &sw, &user));
664fc3eb83cSMatthew G. Knepley     PetscCall(DMGetGlobalVector(dm, &fhat));
665fc3eb83cSMatthew G. Knepley     PetscCall(TestL2ProjectionParticlesToField(dm, sw, fhat, &user));
666fc3eb83cSMatthew G. Knepley     PetscCall(TestL2ProjectionFieldToParticles(dm, sw, fhat, &user));
6679566063dSJacob Faibussowitsch     PetscCall(TestFieldGradientProjection(dm, sw, &user));
668fc3eb83cSMatthew G. Knepley     PetscCall(DMRestoreGlobalVector(dm, &fhat));
66980aecfd0Sjosephpu   } else {
67080aecfd0Sjosephpu     PetscCall(CreateParticles_Shape(dm, &sw, &user));
67180aecfd0Sjosephpu     PetscCall(TestL2Projection_Shape(dm, sw, &user));
67280aecfd0Sjosephpu   }
6739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
6749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
6759566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
676b122ec5aSJacob Faibussowitsch   return 0;
677c4762a1bSJed Brown }
678c4762a1bSJed Brown 
679c4762a1bSJed Brown /*TEST
680c4762a1bSJed Brown 
681832e2b15SMatthew G. Knepley   # Swarm does not handle complex or quad
682832e2b15SMatthew G. Knepley   build:
683832e2b15SMatthew G. Knepley     requires: !complex double
684832e2b15SMatthew G. Knepley 
685c4762a1bSJed Brown   test:
686c4762a1bSJed Brown     suffix: proj_tri_0
687832e2b15SMatthew G. Knepley     requires: triangle
688832e2b15SMatthew 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
689c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
690c4762a1bSJed Brown 
691c4762a1bSJed Brown   test:
692c4762a1bSJed Brown     suffix: proj_tri_2_faces
693832e2b15SMatthew G. Knepley     requires: triangle
694832e2b15SMatthew 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
695c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
696c4762a1bSJed Brown 
697c4762a1bSJed Brown   test:
698c4762a1bSJed Brown     suffix: proj_quad_0
699832e2b15SMatthew G. Knepley     requires: triangle
70030602db0SMatthew 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
701c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
702c4762a1bSJed Brown 
703c4762a1bSJed Brown   test:
704c4762a1bSJed Brown     suffix: proj_quad_2_faces
705832e2b15SMatthew G. Knepley     requires: triangle
70630602db0SMatthew 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
707c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
708*953fc092SRylanor   test:
709*953fc092SRylanor     suffix: proj_quad_hip_serial
710*953fc092SRylanor     requires: triangle hip
711*953fc092SRylanor     nsize: 1
712*953fc092SRylanor     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,8 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -dm_vec_type hip -mat_type aijhipsparse -ptof_ksp_type cg -ptof_pc_type ilu -ptof_ksp_rtol 1e-15 -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
713*953fc092SRylanor     filter: grep -v marker | grep -v atomic | grep -v usage
714*953fc092SRylanor   test:
715*953fc092SRylanor     suffix: proj_quad_hip_parallel
716*953fc092SRylanor     requires: triangle hip
717*953fc092SRylanor     nsize: 2
718*953fc092SRylanor     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,8 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -dm_vec_type hip -mat_type aijhipsparse -ptof_ksp_type cg -ptof_pc_type asm -ptof_sub_pc_type ilu -ptof_ksp_rtol 1e-15 -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
719*953fc092SRylanor     filter: grep -v marker | grep -v atomic | grep -v usage
720c4762a1bSJed Brown   test:
721c4762a1bSJed Brown     suffix: proj_tri_5P
722832e2b15SMatthew G. Knepley     requires: triangle
723832e2b15SMatthew 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
724c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
725c4762a1bSJed Brown 
726c4762a1bSJed Brown   test:
727c4762a1bSJed Brown     suffix: proj_quad_5P
728832e2b15SMatthew G. Knepley     requires: triangle
72930602db0SMatthew 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
730c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
731c4762a1bSJed Brown 
732c4762a1bSJed Brown   test:
733c4762a1bSJed Brown     suffix: proj_tri_mdx
734832e2b15SMatthew G. Knepley     requires: triangle
735832e2b15SMatthew 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
736c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
737c4762a1bSJed Brown 
738c4762a1bSJed Brown   test:
739c4762a1bSJed Brown     suffix: proj_tri_mdx_5P
740832e2b15SMatthew G. Knepley     requires: triangle
741832e2b15SMatthew 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
742c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
743c4762a1bSJed Brown 
744c4762a1bSJed Brown   test:
745c4762a1bSJed Brown     suffix: proj_tri_3d
746832e2b15SMatthew G. Knepley     requires: ctetgen
74730602db0SMatthew 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
748c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
749c4762a1bSJed Brown 
750c4762a1bSJed Brown   test:
751c4762a1bSJed Brown     suffix: proj_tri_3d_2_faces
752832e2b15SMatthew G. Knepley     requires: ctetgen
75330602db0SMatthew 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
754c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
755c4762a1bSJed Brown 
756c4762a1bSJed Brown   test:
757c4762a1bSJed Brown     suffix: proj_tri_3d_5P
758832e2b15SMatthew G. Knepley     requires: ctetgen
75930602db0SMatthew 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
760c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
761c4762a1bSJed Brown 
762c4762a1bSJed Brown   test:
763c4762a1bSJed Brown     suffix: proj_tri_3d_mdx
764832e2b15SMatthew G. Knepley     requires: ctetgen
76530602db0SMatthew 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
766c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
767c4762a1bSJed Brown 
768c4762a1bSJed Brown   test:
769c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_5P
770832e2b15SMatthew G. Knepley     requires: ctetgen
77130602db0SMatthew 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
772c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
773c4762a1bSJed Brown 
774c4762a1bSJed Brown   test:
775c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_2_faces
776832e2b15SMatthew G. Knepley     requires: ctetgen
77730602db0SMatthew 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
778c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
779c4762a1bSJed Brown 
780c4762a1bSJed Brown   test:
781c4762a1bSJed Brown     suffix: proj_tri_3d_mdx_5P_2_faces
782832e2b15SMatthew G. Knepley     requires: ctetgen
78330602db0SMatthew 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
784c4762a1bSJed Brown     filter: grep -v marker | grep -v atomic | grep -v usage
785c4762a1bSJed Brown 
786832e2b15SMatthew G. Knepley   test:
787832e2b15SMatthew G. Knepley     suffix: proj_quad_lsqr_scale
788832e2b15SMatthew G. Knepley     requires: !complex
78930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
790832e2b15SMatthew G. Knepley       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
791832e2b15SMatthew G. Knepley       -particlesPerCell 200 \
792832e2b15SMatthew G. Knepley       -ptof_pc_type lu \
793832e2b15SMatthew G. Knepley       -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none
794832e2b15SMatthew G. Knepley 
795832e2b15SMatthew G. Knepley   test:
796832e2b15SMatthew G. Knepley     suffix: proj_quad_lsqr_prec_scale
797832e2b15SMatthew G. Knepley     requires: !complex
79830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \
799832e2b15SMatthew G. Knepley       -petscspace_degree 2 -petscfe_default_quadrature_order 3 \
800832e2b15SMatthew G. Knepley       -particlesPerCell 200 \
801832e2b15SMatthew G. Knepley       -ptof_pc_type lu \
802fc3eb83cSMatthew 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
80380aecfd0Sjosephpu   test:
80480aecfd0Sjosephpu     suffix: proj_shape_linear_tri_2d
80580aecfd0Sjosephpu     requires: ctetgen triangle
80680aecfd0Sjosephpu     args: -shape -dm_plex_box_faces 2,2\
80780aecfd0Sjosephpu           -dm_view -sw_view\
80880aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
80980aecfd0Sjosephpu           -pc_type lu
81080aecfd0Sjosephpu   test:
81180aecfd0Sjosephpu     suffix: proj_shape_linear_quad_2d
81280aecfd0Sjosephpu     requires: ctetgen triangle
81380aecfd0Sjosephpu     args: -shape -dm_plex_simplex 0 -dm_plex_box_faces 2,2\
81480aecfd0Sjosephpu           -dm_view -sw_view\
81580aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
81680aecfd0Sjosephpu           -pc_type lu
81780aecfd0Sjosephpu   test:
81880aecfd0Sjosephpu     suffix: proj_shape_linear_tri_3d
81980aecfd0Sjosephpu     requires: ctetgen triangle
82080aecfd0Sjosephpu     args: -shape -dm_plex_dim 3 -dm_plex_box_faces 2,2,2\
82180aecfd0Sjosephpu           -dm_view -sw_view\
82280aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
82380aecfd0Sjosephpu           -pc_type lu
82480aecfd0Sjosephpu   test:
82580aecfd0Sjosephpu     suffix: proj_shape_linear_quad_3d
82680aecfd0Sjosephpu     requires: ctetgen triangle
82780aecfd0Sjosephpu     args: -shape -dm_plex_dim 3 -dm_plex_simplex 0\
82880aecfd0Sjosephpu           -dm_plex_box_faces 2,2,2\
82980aecfd0Sjosephpu           -dm_view -sw_view\
83080aecfd0Sjosephpu           -petscspace_degree 1 -petscfe_default_quadrature_order 1\
83180aecfd0Sjosephpu           -pc_type lu
832c4762a1bSJed Brown TEST*/
833