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