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