xref: /petsc/src/dm/impls/swarm/tests/ex4.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1db4653c2SDaniel Finn static char help[] = "Comparing basic symplectic, theta and discrete gradients interators on a simple hamiltonian system (harmonic oscillator) with particles\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>  /* For norm */
5c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */
6c4762a1bSJed Brown #include <petscdmswarm.h>
7c4762a1bSJed Brown #include <petscts.h>
8c4762a1bSJed Brown 
9c4762a1bSJed Brown typedef struct {
10db4653c2SDaniel Finn   PetscInt  dim;                          /* The topological mesh dimension */
11db4653c2SDaniel Finn   PetscBool simplex;                      /* Flag for simplices or tensor cells */
12db4653c2SDaniel Finn   char      filename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */
13c4762a1bSJed Brown   PetscReal omega;                        /* Oscillation frequency omega */
14c4762a1bSJed Brown   PetscInt  particlesPerCell;             /* The number of partices per cell */
15db4653c2SDaniel Finn   PetscInt  numberOfCells;                /* Number of cells in mesh */
16c4762a1bSJed Brown   PetscReal momentTol;                    /* Tolerance for checking moment conservation */
17c4762a1bSJed Brown   PetscBool monitor;                      /* Flag for using the TS monitor */
18c4762a1bSJed Brown   PetscBool error;                        /* Flag for printing the error */
19c4762a1bSJed Brown   PetscInt  ostep;                        /* print the energy at each ostep time steps */
20c4762a1bSJed Brown } AppCtx;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
23c4762a1bSJed Brown {
24c4762a1bSJed Brown   PetscErrorCode ierr;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   PetscFunctionBeginUser;
27db4653c2SDaniel Finn   options->dim              = 2;
28db4653c2SDaniel Finn   options->simplex          = PETSC_TRUE;
29c4762a1bSJed Brown   options->monitor          = PETSC_FALSE;
30c4762a1bSJed Brown   options->error            = PETSC_FALSE;
31c4762a1bSJed Brown   options->particlesPerCell = 1;
32db4653c2SDaniel Finn   options->numberOfCells    = 2;
33c4762a1bSJed Brown   options->momentTol        = 100.0*PETSC_MACHINE_EPSILON;
34c4762a1bSJed Brown   options->omega            = 64.0;
35c4762a1bSJed Brown   options->ostep            = 100;
36c4762a1bSJed Brown 
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcpy(options->filename, ""));
38db4653c2SDaniel Finn 
39c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Harmonic Oscillator Options", "DMPLEX");CHKERRQ(ierr);
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex4.c", options->simplex, &options->simplex, NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex4.c", options->filename, options->filename, sizeof(options->filename), NULL));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-omega", "Oscillator frequency", "ex4.c", options->omega, &options->omega, PETSC_NULL));
48c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   PetscFunctionReturn(0);
51db4653c2SDaniel Finn 
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
55c4762a1bSJed Brown {
56c4762a1bSJed Brown   PetscFunctionBeginUser;
575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
61db4653c2SDaniel Finn 
62c4762a1bSJed Brown   PetscFunctionReturn(0);
63c4762a1bSJed Brown }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw)
66c4762a1bSJed Brown {
67c4762a1bSJed Brown   DM             dm;
68c4762a1bSJed Brown   AppCtx        *user;
69c4762a1bSJed Brown   PetscRandom    rnd;
70c4762a1bSJed Brown   PetscBool      simplex;
71c4762a1bSJed Brown   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
72c4762a1bSJed Brown   PetscInt       dim, d, cStart, cEnd, c, Np, p;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscFunctionBeginUser;
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rnd, -1.0, 1.0));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rnd));
78c4762a1bSJed Brown 
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dmSw, &user));
80c4762a1bSJed Brown   Np   = user->particlesPerCell;
815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(dmSw, &dm));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(dm, &simplex));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
85db4653c2SDaniel Finn   user->numberOfCells = cEnd - cStart;
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
87c4762a1bSJed Brown   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
89c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
90c4762a1bSJed Brown     if (Np == 1) {
915f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
92c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
93c4762a1bSJed Brown     } else {
945f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
95c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
96c4762a1bSJed Brown         const PetscInt n   = c*Np + p;
97c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
98c4762a1bSJed Brown 
99c4762a1bSJed Brown         for (d = 0; d < dim; ++d) {
1005f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscRandomGetValueReal(rnd, &refcoords[d]));
101c4762a1bSJed Brown           sum += refcoords[d];
102c4762a1bSJed Brown         }
103c4762a1bSJed Brown         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
104c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
105c4762a1bSJed Brown       }
106c4762a1bSJed Brown     }
107c4762a1bSJed Brown   }
108db4653c2SDaniel Finn 
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree5(centroid, xi0, v0, J, invJ));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
112c4762a1bSJed Brown   PetscFunctionReturn(0);
113c4762a1bSJed Brown }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
116c4762a1bSJed Brown {
117c4762a1bSJed Brown   DM             dm;
118c4762a1bSJed Brown   AppCtx        *user;
119c4762a1bSJed Brown   PetscReal     *coords;
120c4762a1bSJed Brown   PetscScalar   *initialConditions;
121c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np, p;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBeginUser;
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dmSw, &user));
125c4762a1bSJed Brown   Np   = user->particlesPerCell;
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetCellDM(dmSw, &dm));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(u, &initialConditions));
131c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
132c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
133c4762a1bSJed Brown       const PetscInt n = c*Np + p;
134c4762a1bSJed Brown       initialConditions[n*2+0] = DMPlex_NormD_Internal(dim, &coords[n*dim]);
135c4762a1bSJed Brown       initialConditions[n*2+1] = 0.0;
136c4762a1bSJed Brown     }
137c4762a1bSJed Brown   }
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(u, &initialConditions));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
140c4762a1bSJed Brown   PetscFunctionReturn(0);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
144c4762a1bSJed Brown {
145c4762a1bSJed Brown   PetscInt      *cellid;
146db4653c2SDaniel Finn   PetscInt       dim, cStart, cEnd, c, Np, p;
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   PetscFunctionBeginUser;
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), sw));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*sw, DMSWARM));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(*sw, dim));
153db4653c2SDaniel Finn   Np = user->particlesPerCell;
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(*sw, dm));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(*sw));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*sw));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
162c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
163c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
164c4762a1bSJed Brown       const PetscInt n = c*Np + p;
165c4762a1bSJed Brown 
166c4762a1bSJed Brown       cellid[n] = c;
167c4762a1bSJed Brown     }
168c4762a1bSJed Brown   }
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) *sw, "Particles"));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*sw, NULL, "-sw_view"));
172c4762a1bSJed Brown   PetscFunctionReturn(0);
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
176c4762a1bSJed Brown {
177c4762a1bSJed Brown   AppCtx            *user  = (AppCtx *) ctx;
178c4762a1bSJed Brown   const PetscReal    omega = user->omega;
179c4762a1bSJed Brown   const PetscScalar *u;
180c4762a1bSJed Brown   MPI_Comm           comm;
181c4762a1bSJed Brown   PetscReal          dt;
182c4762a1bSJed Brown   PetscInt           Np, p;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   PetscFunctionBeginUser;
185c4762a1bSJed Brown   if (step%user->ostep == 0) {
1865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetComm((PetscObject) ts, &comm));
1875f80ce2aSJacob Faibussowitsch     if (!step) CHKERRQ(PetscPrintf(comm, "Time     Step Part     Energy Mod Energy\n"));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetTimeStep(ts, &dt));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(U, &u));
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(U, &Np));
191c4762a1bSJed Brown     Np /= 2;
192c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
193c4762a1bSJed Brown       const PetscReal x  = PetscRealPart(u[p*2+0]);
194c4762a1bSJed Brown       const PetscReal v  = PetscRealPart(u[p*2+1]);
195c4762a1bSJed Brown       const PetscReal E  = 0.5*(v*v + PetscSqr(omega)*x*x);
196c4762a1bSJed Brown       const PetscReal mE = 0.5*(v*v + PetscSqr(omega)*x*x - PetscSqr(omega)*dt*x*v);
1975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "%.6lf %4D %4D %10.4lf %10.4lf\n", t, step, p, (double) E, (double) mE));
198c4762a1bSJed Brown     }
1995f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(U, &u));
200c4762a1bSJed Brown   }
201c4762a1bSJed Brown   PetscFunctionReturn(0);
202c4762a1bSJed Brown }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u)
205c4762a1bSJed Brown {
206c4762a1bSJed Brown   DM             dm;
207c4762a1bSJed Brown   AppCtx        *user;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBeginUser;
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dm, &user));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(SetInitialCoordinates(dm));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(SetInitialConditions(dm, u));
214c4762a1bSJed Brown   PetscFunctionReturn(0);
215c4762a1bSJed Brown }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
218c4762a1bSJed Brown {
219c4762a1bSJed Brown   MPI_Comm           comm;
220c4762a1bSJed Brown   DM                 sdm;
221c4762a1bSJed Brown   AppCtx            *user;
222c4762a1bSJed Brown   const PetscScalar *u, *coords;
223c4762a1bSJed Brown   PetscScalar       *e;
224c4762a1bSJed Brown   PetscReal          t, omega;
225c4762a1bSJed Brown   PetscInt           dim, Np, p;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   PetscFunctionBeginUser;
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) ts, &comm));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &sdm));
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(sdm, &user));
231c4762a1bSJed Brown   omega = user->omega;
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(sdm, &dim));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts, &t));
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(E, &e));
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U, &u));
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(U, &Np));
237c4762a1bSJed Brown   Np  /= 2;
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
239c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
240c4762a1bSJed Brown     const PetscReal x  = PetscRealPart(u[p*2+0]);
241c4762a1bSJed Brown     const PetscReal v  = PetscRealPart(u[p*2+1]);
242c4762a1bSJed Brown     const PetscReal x0 = DMPlex_NormD_Internal(dim, &coords[p*dim]);
243c4762a1bSJed Brown     const PetscReal ex =  x0*PetscCosReal(omega*t);
244c4762a1bSJed Brown     const PetscReal ev = -x0*omega*PetscSinReal(omega*t);
245c4762a1bSJed Brown 
2465f80ce2aSJacob Faibussowitsch     if (user->error) CHKERRQ(PetscPrintf(comm, "p%D error [%.2g %.2g] sol [%.6lf %.6lf] exact [%.6lf %.6lf] energy/exact energy %g / %g\n", p, (double) PetscAbsReal(x-ex), (double) PetscAbsReal(v-ev), (double) x, (double) v, (double) ex, (double) ev, 0.5*(v*v + PetscSqr(omega)*x*x), (double) 0.5*PetscSqr(omega*x0)));
247c4762a1bSJed Brown     e[p*2+0] = x - ex;
248c4762a1bSJed Brown     e[p*2+1] = v - ev;
249c4762a1bSJed Brown   }
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
25140be0ff1SMatthew G. Knepley 
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U, &u));
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(E, &e));
254c4762a1bSJed Brown   PetscFunctionReturn(0);
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
25740be0ff1SMatthew G. Knepley /*---------------------Create particle RHS Functions--------------------------*/
25840be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
25940be0ff1SMatthew G. Knepley {
26040be0ff1SMatthew G. Knepley   const PetscScalar *v;
26140be0ff1SMatthew G. Knepley   PetscScalar       *xres;
26240be0ff1SMatthew G. Knepley   PetscInt          Np, p;
26340be0ff1SMatthew G. Knepley 
26440be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(Xres, &xres));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(V, &v));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(Xres, &Np));
26840be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
26940be0ff1SMatthew G. Knepley      xres[p] = v[p];
27040be0ff1SMatthew G. Knepley   }
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(V, &v));
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(Xres, &xres));
27340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
27440be0ff1SMatthew G. Knepley }
27540be0ff1SMatthew G. Knepley 
27640be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
27740be0ff1SMatthew G. Knepley {
27840be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *)ctx;
27940be0ff1SMatthew G. Knepley   const PetscScalar *x;
28040be0ff1SMatthew G. Knepley   PetscScalar       *vres;
28140be0ff1SMatthew G. Knepley   PetscInt          Np, p;
28240be0ff1SMatthew G. Knepley 
28340be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(Vres, &vres));
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X, &x));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(Vres, &Np));
28740be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
28840be0ff1SMatthew G. Knepley     vres[p] = -PetscSqr(user->omega)*x[p];
28940be0ff1SMatthew G. Knepley   }
2905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X, &x));
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(Vres, &vres));
29240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
29340be0ff1SMatthew G. Knepley }
29440be0ff1SMatthew G. Knepley 
29540be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/
29640be0ff1SMatthew G. Knepley 
29740be0ff1SMatthew G. Knepley /*--------------------Define RHSFunction, RHSJacobian (Theta)-----------------*/
29840be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
29940be0ff1SMatthew G. Knepley {
30040be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
30140be0ff1SMatthew G. Knepley   DM                dm;
30240be0ff1SMatthew G. Knepley   const PetscScalar *u;
30340be0ff1SMatthew G. Knepley   PetscScalar       *g;
30440be0ff1SMatthew G. Knepley   PetscInt          Np, p;
30540be0ff1SMatthew G. Knepley 
30640be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
3075f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(G, &g));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U, &u));
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(U, &Np));
31140be0ff1SMatthew G. Knepley   Np  /= 2;
31240be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
31340be0ff1SMatthew G. Knepley     g[p*2+0] = u[p*2+1];
31440be0ff1SMatthew G. Knepley     g[p*2+1] = -PetscSqr(user->omega)*u[p*2+0];
31540be0ff1SMatthew G. Knepley   }
3165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U, &u));
3175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(G, &g));
31840be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
31940be0ff1SMatthew G. Knepley }
32040be0ff1SMatthew G. Knepley 
32140be0ff1SMatthew G. Knepley /*Ji = dFi/dxj
32240be0ff1SMatthew G. Knepley J= (0    1)
32340be0ff1SMatthew G. Knepley    (-w^2 0)
32440be0ff1SMatthew G. Knepley */
32540be0ff1SMatthew G. Knepley static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U , Mat J, Mat P, void *ctx)
32640be0ff1SMatthew G. Knepley {
32740be0ff1SMatthew G. Knepley   AppCtx             *user = (AppCtx *) ctx;
328db4653c2SDaniel Finn   PetscInt           Np;
32940be0ff1SMatthew G. Knepley   PetscInt           i, m, n;
33040be0ff1SMatthew G. Knepley   const PetscScalar *u;
33140be0ff1SMatthew G. Knepley   PetscScalar        vals[4] = {0., 1., -PetscSqr(user->omega), 0.};
33240be0ff1SMatthew G. Knepley 
333362febeeSStefano Zampini   PetscFunctionBeginUser;
3345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U, &u));
3355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(U, &Np));
336db4653c2SDaniel Finn   Np /= 2;
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(J, &m, &n));
33840be0ff1SMatthew G. Knepley   for (i = 0; i < Np; ++i) {
33940be0ff1SMatthew G. Knepley     const PetscInt rows[2] = {2*i, 2*i+1};
3405f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
34140be0ff1SMatthew G. Knepley   }
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
3435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
34440be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
34540be0ff1SMatthew G. Knepley 
34640be0ff1SMatthew G. Knepley }
347db4653c2SDaniel Finn 
34840be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/
34940be0ff1SMatthew G. Knepley 
35040be0ff1SMatthew G. Knepley /*----------------Define S, F, G Functions (Discrete Gradients)---------------*/
35140be0ff1SMatthew G. Knepley /*
35240be0ff1SMatthew G. Knepley   u_t = S * gradF
35340be0ff1SMatthew G. Knepley     --or--
35440be0ff1SMatthew G. Knepley   u_t = S * G
35540be0ff1SMatthew G. Knepley 
35640be0ff1SMatthew G. Knepley   + Sfunc - constructor for the S matrix from the formulation
35740be0ff1SMatthew G. Knepley   . Ffunc - functional F from the formulation
35840be0ff1SMatthew G. Knepley   - Gfunc - constructor for the gradient of F from the formulation
35940be0ff1SMatthew G. Knepley */
36040be0ff1SMatthew G. Knepley 
36140be0ff1SMatthew G. Knepley PetscErrorCode Sfunc(TS ts, PetscReal t, Vec U, Mat S, void *ctx)
36240be0ff1SMatthew G. Knepley {
363db4653c2SDaniel Finn   PetscInt           Np;
36440be0ff1SMatthew G. Knepley   PetscInt           i, m, n;
36540be0ff1SMatthew G. Knepley   const PetscScalar *u;
36640be0ff1SMatthew G. Knepley   PetscScalar        vals[4] = {0., 1., -1, 0.};
36740be0ff1SMatthew G. Knepley 
36840be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
3695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U, &u));
3705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(U, &Np));
371db4653c2SDaniel Finn   Np /= 2;
3725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(S, &m, &n));
37340be0ff1SMatthew G. Knepley   for (i = 0; i < Np; ++i) {
37440be0ff1SMatthew G. Knepley     const PetscInt rows[2] = {2*i, 2*i+1};
3755f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES));
37640be0ff1SMatthew G. Knepley   }
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY));
3785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY));
37940be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
380db4653c2SDaniel Finn 
38140be0ff1SMatthew G. Knepley }
38240be0ff1SMatthew G. Knepley 
38340be0ff1SMatthew G. Knepley PetscErrorCode Ffunc(TS ts, PetscReal t, Vec U, PetscScalar *F, void *ctx)
38440be0ff1SMatthew G. Knepley {
38540be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
38640be0ff1SMatthew G. Knepley   DM                 dm;
38740be0ff1SMatthew G. Knepley   const PetscScalar *u;
388db4653c2SDaniel Finn   PetscInt           Np;
38940be0ff1SMatthew G. Knepley   PetscInt           p;
39040be0ff1SMatthew G. Knepley 
39140be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
3925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
393db4653c2SDaniel Finn 
394db4653c2SDaniel Finn   /* Define F */
3955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U, &u));
3965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(U, &Np));
397db4653c2SDaniel Finn   Np /= 2;
39840be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
399db4653c2SDaniel Finn     *F += 0.5*PetscSqr(user->omega)*PetscSqr(u[p*2+0]) + 0.5*PetscSqr(u[p*2+1]);
40040be0ff1SMatthew G. Knepley   }
4015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U, &u));
40240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
40340be0ff1SMatthew G. Knepley }
40440be0ff1SMatthew G. Knepley 
40540be0ff1SMatthew G. Knepley PetscErrorCode gradFfunc(TS ts, PetscReal t, Vec U, Vec gradF, void *ctx)
40640be0ff1SMatthew G. Knepley {
40740be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
40840be0ff1SMatthew G. Knepley   DM                dm;
40940be0ff1SMatthew G. Knepley   const PetscScalar *u;
41040be0ff1SMatthew G. Knepley   PetscScalar       *g;
411db4653c2SDaniel Finn   PetscInt          Np;
41240be0ff1SMatthew G. Knepley   PetscInt          p;
41340be0ff1SMatthew G. Knepley 
41440be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
4155f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
4165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U, &u));
4175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(U, &Np));
418db4653c2SDaniel Finn   Np /= 2;
419db4653c2SDaniel Finn   /*Define gradF*/
4205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(gradF, &g));
42140be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
42240be0ff1SMatthew G. Knepley     g[p*2+0] = PetscSqr(user->omega)*u[p*2+0]; /*dF/dx*/
42340be0ff1SMatthew G. Knepley     g[p*2+1] = u[p*2+1]; /*dF/dv*/
42440be0ff1SMatthew G. Knepley   }
4255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U, &u));
4265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(gradF, &g));
42740be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
42840be0ff1SMatthew G. Knepley }
429db4653c2SDaniel Finn 
43040be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/
43140be0ff1SMatthew G. Knepley 
432c4762a1bSJed Brown int main(int argc,char **argv)
433c4762a1bSJed Brown {
434c4762a1bSJed Brown   TS             ts;     /* nonlinear solver                 */
435c4762a1bSJed Brown   DM             dm, sw; /* Mesh and particle managers       */
436c4762a1bSJed Brown   Vec            u;      /* swarm vector                     */
43740be0ff1SMatthew G. Knepley   PetscInt       n;      /* swarm vector size                */
438c4762a1bSJed Brown   IS             is1, is2;
439c4762a1bSJed Brown   MPI_Comm       comm;
44040be0ff1SMatthew G. Knepley   Mat            J;      /* Jacobian matrix                  */
441c4762a1bSJed Brown   AppCtx         user;
442c4762a1bSJed Brown 
44340be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44440be0ff1SMatthew G. Knepley      Initialize program
44540be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
446*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
447c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
4485f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(comm, &user));
449c4762a1bSJed Brown 
45040be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45140be0ff1SMatthew G. Knepley      Create Particle-Mesh
45240be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4535f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(comm, &dm, &user));
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateParticles(dm, &sw, &user));
4555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(sw, &user));
456c4762a1bSJed Brown 
45740be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45840be0ff1SMatthew G. Knepley      Setup timestepping solver
45940be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(comm, &ts));
4615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
46240be0ff1SMatthew G. Knepley 
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts, sw));
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts, 0.1));
4655f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts, 0.00001));
4665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts, 100));
4675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
4685f80ce2aSJacob Faibussowitsch   if (user.monitor) CHKERRQ(TSMonitorSet(ts, Monitor, &user, NULL));
469c4762a1bSJed Brown 
47040be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47140be0ff1SMatthew G. Knepley      Prepare to solve
47240be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u));
4745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(u, &n));
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
4765f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetComputeInitialCondition(ts, InitializeSolve));
4775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetComputeExactError(ts, ComputeError));
47840be0ff1SMatthew G. Knepley 
47940be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48040be0ff1SMatthew G. Knepley      Define function F(U, Udot , x , t) = G(U, x, t)
48140be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48240be0ff1SMatthew G. Knepley 
48340be0ff1SMatthew G. Knepley   /* - - - - - - - Basic Symplectic - - - - - - - - - - - - - - - - - - - - - -*/
4845f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(comm, n/2, 0, 2, &is1));
4855f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(comm, n/2, 1, 2, &is2));
4865f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts, "position", is1));
4875f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetIS(ts, "momentum", is2));
4885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
4895f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is2));
4905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunction1, &user));
4915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunction2, &user));
492c4762a1bSJed Brown 
49340be0ff1SMatthew G. Knepley   /* - - - - - - - Theta (Implicit Midpoint) - - - - - - - - - - - - - - - - -*/
494db4653c2SDaniel Finn 
4955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
49640be0ff1SMatthew G. Knepley 
4975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&J));
4985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(J));
5005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(J));
5015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
5025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
5035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,J,J,RHSJacobian,&user));
50440be0ff1SMatthew G. Knepley 
50540be0ff1SMatthew G. Knepley   /* - - - - - - - Discrete Gradients - - - - - - - - - - - - - - - - - - - - */
5065f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDiscGradSetFormulation(ts, Sfunc, Ffunc, gradFfunc, &user));
50740be0ff1SMatthew G. Knepley 
50840be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50940be0ff1SMatthew G. Knepley      Solve
51040be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
5115f80ce2aSJacob Faibussowitsch   CHKERRQ(TSComputeInitialCondition(ts, u));
5125f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts, u));
513c4762a1bSJed Brown 
51440be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51540be0ff1SMatthew G. Knepley      Clean up workspace
51640be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
5175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
5185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u));
5195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
5205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&sw));
5215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
522*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
523*b122ec5aSJacob Faibussowitsch   return 0;
524c4762a1bSJed Brown }
525c4762a1bSJed Brown 
526c4762a1bSJed Brown /*TEST
527c4762a1bSJed Brown 
528c4762a1bSJed Brown    build:
529c4762a1bSJed Brown      requires: triangle !single !complex
530c4762a1bSJed Brown    test:
531c4762a1bSJed Brown      suffix: 1
532db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 1 -ts_convergence_estimate -convest_num_refine 2 -dm_view  -monitor -output_step 50 -error
53340be0ff1SMatthew G. Knepley 
534c4762a1bSJed Brown    test:
535c4762a1bSJed Brown      suffix: 2
536db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 2 -ts_convergence_estimate -convest_num_refine 2 -dm_view  -monitor -output_step 50 -error
53740be0ff1SMatthew G. Knepley 
538c4762a1bSJed Brown    test:
539c4762a1bSJed Brown      suffix: 3
540db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -monitor -output_step 50 -error
54140be0ff1SMatthew G. Knepley 
542c4762a1bSJed Brown    test:
543c4762a1bSJed Brown      suffix: 4
544db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 4 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view  -monitor -output_step 50 -error
54540be0ff1SMatthew G. Knepley 
54640be0ff1SMatthew G. Knepley    test:
54740be0ff1SMatthew G. Knepley      suffix: 5
548db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type theta -ts_theta_theta 0.5 -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2 -dm_view
54940be0ff1SMatthew G. Knepley 
55040be0ff1SMatthew G. Knepley    test:
55140be0ff1SMatthew G. Knepley      suffix: 6
552db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type discgrad -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2 -dm_view
553db4653c2SDaniel Finn 
554db4653c2SDaniel Finn    test:
555db4653c2SDaniel Finn      suffix: 7
556db4653c2SDaniel Finn      args: -dm_plex_box_faces 1,1 -ts_type discgrad -ts_discgrad_gonzalez -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2 -dm_view
557c4762a1bSJed Brown 
558c4762a1bSJed Brown TEST*/
559