xref: /petsc/src/ts/tutorials/hamiltonian/ex2.c (revision d0609ced746bc51b019815ca91d747429db24893)
1b80bf5b1SJoe Pusztay static char help[] = "Two stream instability from Birdsal and Langdon with DMSwarm and TS basic symplectic integrators\n";
2b80bf5b1SJoe Pusztay 
3b80bf5b1SJoe Pusztay #include <petscdmplex.h>
4b80bf5b1SJoe Pusztay #include <petscfe.h>
5b80bf5b1SJoe Pusztay #include <petscdmswarm.h>
6b80bf5b1SJoe Pusztay #include <petscds.h>
7b80bf5b1SJoe Pusztay #include <petscksp.h>
8b80bf5b1SJoe Pusztay #include <petsc/private/petscfeimpl.h>
9b80bf5b1SJoe Pusztay #include <petsc/private/tsimpl.h>
10b80bf5b1SJoe Pusztay #include <petscts.h>
11b80bf5b1SJoe Pusztay #include <petscmath.h>
12b80bf5b1SJoe Pusztay 
13b80bf5b1SJoe Pusztay typedef struct {
14b80bf5b1SJoe Pusztay   PetscInt       dim;                              /* The topological mesh dimension */
15b80bf5b1SJoe Pusztay   PetscBool      simplex;                          /* Flag for simplices or tensor cells */
16b80bf5b1SJoe Pusztay   PetscBool      bdm;                              /* Flag for mixed form poisson */
17b80bf5b1SJoe Pusztay   PetscBool      monitor;                          /* Flag for use of the TS monitor */
18b80bf5b1SJoe Pusztay   PetscBool      uniform;                          /* Flag to uniformly space particles in x */
19b80bf5b1SJoe Pusztay   char           meshFilename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */
20b80bf5b1SJoe Pusztay   PetscReal      sigma;                            /* Linear charge per box length */
21b80bf5b1SJoe Pusztay   PetscReal      timeScale;                        /* Nondimensionalizing time scaling */
22b80bf5b1SJoe Pusztay   PetscInt       particlesPerCell;                 /* The number of partices per cell */
23b80bf5b1SJoe Pusztay   PetscReal      particleRelDx;                    /* Relative particle position perturbation compared to average cell diameter h */
24b80bf5b1SJoe Pusztay   PetscInt       k;                                /* Mode number for test function */
25b80bf5b1SJoe Pusztay   PetscReal      momentTol;                        /* Tolerance for checking moment conservation */
26b80bf5b1SJoe Pusztay   SNES           snes;                             /* SNES object */
27b80bf5b1SJoe Pusztay   PetscInt       steps;                            /* TS iterations */
28b80bf5b1SJoe Pusztay   PetscReal      stepSize;                         /* Time stepper step size */
29b80bf5b1SJoe Pusztay   PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
30b80bf5b1SJoe Pusztay } AppCtx;
31b80bf5b1SJoe Pusztay 
32b80bf5b1SJoe Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33b80bf5b1SJoe Pusztay {
34b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
35b80bf5b1SJoe Pusztay   options->dim              = 2;
36b80bf5b1SJoe Pusztay   options->simplex          = PETSC_TRUE;
37b80bf5b1SJoe Pusztay   options->monitor          = PETSC_TRUE;
38b80bf5b1SJoe Pusztay   options->particlesPerCell = 1;
39b80bf5b1SJoe Pusztay   options->k                = 1;
40b80bf5b1SJoe Pusztay   options->particleRelDx    = 1.e-20;
41b80bf5b1SJoe Pusztay   options->momentTol        = 100.*PETSC_MACHINE_EPSILON;
42b80bf5b1SJoe Pusztay   options->sigma            = 1.;
43b80bf5b1SJoe Pusztay   options->timeScale        = 1.0e-6;
44b80bf5b1SJoe Pusztay   options->uniform          = PETSC_FALSE;
45b80bf5b1SJoe Pusztay   options->steps            = 1;
46b80bf5b1SJoe Pusztay   options->stepSize         = 0.01;
47b80bf5b1SJoe Pusztay   options->bdm              = PETSC_FALSE;
48b80bf5b1SJoe Pusztay 
49*d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Two Stream options", "DMPLEX");
509566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(options->meshFilename, ""));
519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dim", "The topological mesh dimension", "ex2.c", options->dim, &options->dim, NULL));
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-steps", "TS steps to take", "ex2.c", options->steps, &options->steps, NULL));
539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitor", "To use the TS monitor or not", "ex2.c", options->monitor, &options->monitor, NULL));
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex2.c", options->simplex, &options->simplex, NULL));
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-uniform", "Uniform particle spacing", "ex2.c", options->uniform, &options->uniform, NULL));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-bdm", "Use H1 instead of C0", "ex2.c", options->bdm, &options->bdm, NULL));
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex2.c", options->meshFilename, options->meshFilename, PETSC_MAX_PATH_LEN, NULL));
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex5.c", options->k, &options->k, NULL));
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-sigma","parameter","<1>",options->sigma,&options->sigma,PETSC_NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-stepSize","parameter","<1e-2>",options->stepSize,&options->stepSize,PETSC_NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-timeScale","parameter","<1>",options->timeScale,&options->timeScale,PETSC_NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL));
64*d0609cedSBarry Smith   PetscOptionsEnd();
65b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
66b80bf5b1SJoe Pusztay }
67b80bf5b1SJoe Pusztay 
68b80bf5b1SJoe Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
69b80bf5b1SJoe Pusztay {
70b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
719566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
729566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
739566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
749566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
75b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
76b80bf5b1SJoe Pusztay }
77b80bf5b1SJoe Pusztay 
78b80bf5b1SJoe Pusztay static void laplacian_f1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
79b80bf5b1SJoe Pusztay                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
80b80bf5b1SJoe Pusztay                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
81b80bf5b1SJoe Pusztay                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
82b80bf5b1SJoe Pusztay {
83b80bf5b1SJoe Pusztay   PetscInt d;
84b80bf5b1SJoe Pusztay   for (d = 0; d < dim; ++d) {f1[d] = u_x[d];}
85b80bf5b1SJoe Pusztay }
86b80bf5b1SJoe Pusztay 
87b80bf5b1SJoe Pusztay static void laplacian(PetscInt dim, PetscInt Nf, PetscInt NfAux,
88b80bf5b1SJoe Pusztay                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
89b80bf5b1SJoe Pusztay                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90b80bf5b1SJoe Pusztay                       PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
91b80bf5b1SJoe Pusztay {
92b80bf5b1SJoe Pusztay   PetscInt d;
93b80bf5b1SJoe Pusztay   for (d = 0; d < dim; ++d) {g3[d*dim+d] = 1.0;}
94b80bf5b1SJoe Pusztay }
95b80bf5b1SJoe Pusztay 
96b80bf5b1SJoe Pusztay static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
97b80bf5b1SJoe Pusztay {
98b80bf5b1SJoe Pusztay   PetscFE        fe;
99b80bf5b1SJoe Pusztay   PetscDS        ds;
100b80bf5b1SJoe Pusztay   DMPolytopeType ct;
101b80bf5b1SJoe Pusztay   PetscBool      simplex;
102b80bf5b1SJoe Pusztay   PetscInt       dim, cStart;
103b80bf5b1SJoe Pusztay 
104b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
1059566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
108b80bf5b1SJoe Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1099566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, NULL, -1, &fe));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, "potential"));
1119566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
1129566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1139566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1149566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1159566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, NULL, laplacian_f1));
1169566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, laplacian));
117b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
118b80bf5b1SJoe Pusztay }
119b80bf5b1SJoe Pusztay 
120b80bf5b1SJoe Pusztay /*
121b80bf5b1SJoe Pusztay   Initialize particle coordinates uniformly and with opposing velocities
122b80bf5b1SJoe Pusztay */
123b80bf5b1SJoe Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
124b80bf5b1SJoe Pusztay {
125b80bf5b1SJoe Pusztay   PetscRandom    rnd, rndp;
126b80bf5b1SJoe Pusztay   PetscReal      interval = user->particleRelDx;
127b80bf5b1SJoe Pusztay   PetscScalar    value, *vals;
128b80bf5b1SJoe Pusztay   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *initialConditions, normalized_vel;
129b80bf5b1SJoe Pusztay   PetscInt      *cellid, cStart;
130b80bf5b1SJoe Pusztay   PetscInt       Ncell, Np = user->particlesPerCell, p, c, dim, d;
131b80bf5b1SJoe Pusztay 
132b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
1339566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1349566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
1359566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1369566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
1379566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
1389566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, 0.0, 1.0));
1399566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
1409566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rndp));
1419566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rndp, -interval, interval));
1429566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rndp));
1439566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1449566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
1459566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1469566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL));
1479566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
1489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &Ncell));
1499566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0));
1509566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
1519566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
1529566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
1539566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **) &vals));
1549566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, "kinematics", NULL, NULL, (void **) &initialConditions));
1559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
156b80bf5b1SJoe Pusztay   for (c = cStart; c < Ncell; c++) {
157b80bf5b1SJoe Pusztay     if (Np == 1) {
1589566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
159b80bf5b1SJoe Pusztay       cellid[c] = c;
160b80bf5b1SJoe Pusztay       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
161b80bf5b1SJoe Pusztay     } else {
162b80bf5b1SJoe Pusztay       for (d = 0; d < dim; ++d) xi0[d] = -1.0;
1639566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
164b80bf5b1SJoe Pusztay       for (p = 0; p < Np; ++p) {
165b80bf5b1SJoe Pusztay         const PetscInt n   = c*Np + p;
166b80bf5b1SJoe Pusztay         PetscReal      refcoords[3], spacing;
167b80bf5b1SJoe Pusztay 
168b80bf5b1SJoe Pusztay         cellid[n] = c;
169b80bf5b1SJoe Pusztay         if (user->uniform) {
170b80bf5b1SJoe Pusztay           spacing = 2./Np;
1719566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValue(rnd, &value));
172b80bf5b1SJoe Pusztay           for (d=0; d<dim; ++d) refcoords[d] = d == 0 ? -1. + spacing/2. + p*spacing + value/100. : 0.;
173b80bf5b1SJoe Pusztay         }
174b80bf5b1SJoe Pusztay         else{
1759566063dSJacob Faibussowitsch           for (d = 0; d < dim; ++d) {PetscCall(PetscRandomGetValue(rnd, &value)); refcoords[d] = d == 0 ? PetscRealPart(value) : 0. ;}
176b80bf5b1SJoe Pusztay         }
177b80bf5b1SJoe Pusztay         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
178b80bf5b1SJoe Pusztay         /* constant particle weights */
179b80bf5b1SJoe Pusztay         for (d = 0; d < dim; ++d) vals[n] = user->sigma/Np;
180b80bf5b1SJoe Pusztay       }
181b80bf5b1SJoe Pusztay     }
182b80bf5b1SJoe Pusztay   }
1839566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
184b80bf5b1SJoe Pusztay   normalized_vel = 1.;
185b80bf5b1SJoe Pusztay   for (c = 0; c < Ncell; ++c) {
186b80bf5b1SJoe Pusztay     for (p = 0; p < Np; ++p) {
187b80bf5b1SJoe Pusztay       if (p%2 == 0) {
188b80bf5b1SJoe Pusztay         for (d = 0; d < dim; ++d) initialConditions[(c*Np + p)*dim + d] = d == 0 ? normalized_vel : 0.;
189b80bf5b1SJoe Pusztay       }
190b80bf5b1SJoe Pusztay       else {
191b80bf5b1SJoe Pusztay         for (d = 0; d < dim; ++d) initialConditions[(c*Np + p)*dim + d] = d == 0 ? -(normalized_vel) : 0.;
192b80bf5b1SJoe Pusztay       }
193b80bf5b1SJoe Pusztay     }
194b80bf5b1SJoe Pusztay   }
1959566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
1969566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
1979566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **) &vals));
1989566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, "kinematics", NULL, NULL, (void **) &initialConditions));
1999566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
2009566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rndp));
2019566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
2029566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
2039566063dSJacob Faibussowitsch   PetscCall(DMLocalizeCoordinates(*sw));
204b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
205b80bf5b1SJoe Pusztay }
206b80bf5b1SJoe Pusztay 
207b80bf5b1SJoe Pusztay /* Solve for particle position updates */
208b80bf5b1SJoe Pusztay static PetscErrorCode RHSFunction1(TS ts,PetscReal t,Vec V,Vec Posres,void *ctx)
209b80bf5b1SJoe Pusztay {
210b80bf5b1SJoe Pusztay   const PetscScalar *v;
211b80bf5b1SJoe Pusztay   PetscScalar       *posres;
212b80bf5b1SJoe Pusztay   PetscInt          Np, p, dim, d;
213b80bf5b1SJoe Pusztay   DM                dm;
214b80bf5b1SJoe Pusztay 
215b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
2169566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Posres, &Np));
2179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Posres,&posres));
2189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V,&v));
2199566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2209566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
221b80bf5b1SJoe Pusztay   Np  /= dim;
222b80bf5b1SJoe Pusztay   for (p = 0; p < Np; ++p) {
223b80bf5b1SJoe Pusztay      for (d = 0; d < dim; ++d) {
224b80bf5b1SJoe Pusztay        posres[p*dim+d] = v[p*dim+d];
225b80bf5b1SJoe Pusztay      }
226b80bf5b1SJoe Pusztay   }
2279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V,&v));
2289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Posres,&posres));
229b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
230b80bf5b1SJoe Pusztay 
231b80bf5b1SJoe Pusztay }
232b80bf5b1SJoe Pusztay 
233b80bf5b1SJoe Pusztay /*
234b80bf5b1SJoe Pusztay   Solve for the gradient of the electric field and apply force to particles.
235b80bf5b1SJoe Pusztay  */
236b80bf5b1SJoe Pusztay static PetscErrorCode RHSFunction2(TS ts,PetscReal t,Vec X,Vec Vres,void *ctx)
237b80bf5b1SJoe Pusztay {
238b80bf5b1SJoe Pusztay  AppCtx             *user = (AppCtx *) ctx;
239b80bf5b1SJoe Pusztay   DM                 dm, plex;
240b80bf5b1SJoe Pusztay   PetscDS            prob;
241b80bf5b1SJoe Pusztay   PetscFE            fe;
242b80bf5b1SJoe Pusztay   Mat                M_p;
243b80bf5b1SJoe Pusztay   Vec                phi, locPhi, rho, f;
244b80bf5b1SJoe Pusztay   const PetscScalar *x;
245b80bf5b1SJoe Pusztay   PetscScalar       *vres;
246b80bf5b1SJoe Pusztay   PetscReal         *coords, phi_0;
247b80bf5b1SJoe Pusztay   PetscInt           dim, d, cStart, cEnd, cell, cdim;
248b80bf5b1SJoe Pusztay   PetscReal          m_e = 9.11e-31, q_e = 1.60e-19, epsi_0 = 8.85e-12;
249b80bf5b1SJoe Pusztay 
250b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
251b80bf5b1SJoe Pusztay   PetscObjectSetName((PetscObject) X, "rhsf2 position");
252b80bf5b1SJoe Pusztay   VecViewFromOptions(X, NULL, "-rhsf2_x_view");
2539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vres,&vres));
2559566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2569566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2579566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(user->snes, &plex));
2589566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(plex, &cdim));
2599566063dSJacob Faibussowitsch   PetscCall(DMGetDS(plex, &prob));
2609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe));
2619566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(plex, &phi));
2629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(plex, &locPhi));
2639566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, plex, &M_p));
2649566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
2659566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(plex, &rho));
2669566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, "w_q", &f));
2679566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) f, "weights vector"));
2689566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
2699566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(M_p, f, rho));
2709566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, "w_q", &f));
2719566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) rho, "rho"));
2729566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
273b80bf5b1SJoe Pusztay   /* Take nullspace out of rhs */
274b80bf5b1SJoe Pusztay   {
275b80bf5b1SJoe Pusztay     PetscScalar sum;
276b80bf5b1SJoe Pusztay     PetscInt    n;
277b80bf5b1SJoe Pusztay     phi_0 = (user->sigma*user->sigma*user->sigma)*(user->timeScale*user->timeScale)/(m_e*q_e*epsi_0);
278b80bf5b1SJoe Pusztay 
2799566063dSJacob Faibussowitsch     PetscCall(VecGetSize(rho, &n));
2809566063dSJacob Faibussowitsch     PetscCall(VecSum(rho, &sum));
2819566063dSJacob Faibussowitsch     PetscCall(VecShift(rho, -sum/n));
282b80bf5b1SJoe Pusztay 
2839566063dSJacob Faibussowitsch     PetscCall(VecSum(rho, &sum));
2843c633725SBarry Smith     PetscCheck(PetscAbsScalar(sum) <= 1.0e-10,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Charge should have no DC component %g", sum);
2859566063dSJacob Faibussowitsch     PetscCall(VecScale(rho, phi_0));
286b80bf5b1SJoe Pusztay   }
2879566063dSJacob Faibussowitsch   PetscCall(VecSet(phi, 0.0));
2889566063dSJacob Faibussowitsch   PetscCall(SNESSolve(user->snes, rho, phi));
2899566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
2909566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(plex, &rho));
2919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M_p));
2929566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(plex, phi, INSERT_VALUES, locPhi));
2939566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(plex, phi, INSERT_VALUES, locPhi));
2949566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dm));
2959566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
2969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
297b80bf5b1SJoe Pusztay   for (cell = cStart; cell < cEnd; ++cell) {
298b80bf5b1SJoe Pusztay     PetscTabulation tab;
299b80bf5b1SJoe Pusztay     PetscReal    v[3], J[9], invJ[9], detJ;
300b80bf5b1SJoe Pusztay     PetscScalar *ph       = PETSC_NULL;
301b80bf5b1SJoe Pusztay     PetscReal   *pcoord   = PETSC_NULL;
302b80bf5b1SJoe Pusztay     PetscReal   *refcoord = PETSC_NULL;
303b80bf5b1SJoe Pusztay     PetscInt    *points   = PETSC_NULL, Ncp, cp;
304b80bf5b1SJoe Pusztay     PetscScalar  gradPhi[3];
305b80bf5b1SJoe Pusztay 
3069566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, NULL, v, J, invJ, &detJ));
3079566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dm, cell, &Ncp, &points));
3089566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, Ncp*cdim, MPIU_REAL, &pcoord));
3099566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, Ncp*cdim, MPIU_REAL, &refcoord));
310b80bf5b1SJoe Pusztay     for (cp = 0; cp < Ncp; ++cp) {
311b80bf5b1SJoe Pusztay       for (d = 0; d < cdim; ++d) {
312b80bf5b1SJoe Pusztay         pcoord[cp*cdim+d] = coords[points[cp]*cdim+d];
313b80bf5b1SJoe Pusztay       }
314b80bf5b1SJoe Pusztay     }
3159566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesToReference(plex, cell, Ncp, pcoord, refcoord));
3169566063dSJacob Faibussowitsch     PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
3179566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, NULL, locPhi, cell, NULL, &ph));
318b80bf5b1SJoe Pusztay     for (cp = 0; cp < Ncp; ++cp) {
319b80bf5b1SJoe Pusztay       const PetscInt p = points[cp];
320b80bf5b1SJoe Pusztay       gradPhi[0] = 0.0;
321b80bf5b1SJoe Pusztay       gradPhi[1] = 0.0;
322b80bf5b1SJoe Pusztay       gradPhi[2] = 0.0;
323b80bf5b1SJoe Pusztay       const PetscReal *basisDer = tab->T[1];
324b80bf5b1SJoe Pusztay 
3259566063dSJacob Faibussowitsch       PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, ph, cdim, invJ, NULL, cp, gradPhi));
326b80bf5b1SJoe Pusztay       for (d = 0; d < cdim; ++d) {
327b80bf5b1SJoe Pusztay         vres[p*cdim+d] = d == 0 ? gradPhi[d] : 0.;
328b80bf5b1SJoe Pusztay       }
329b80bf5b1SJoe Pusztay     }
3309566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(plex, NULL, locPhi, cell, NULL, &ph));
3319566063dSJacob Faibussowitsch     PetscCall(PetscTabulationDestroy(&tab));
3329566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Ncp*cdim, MPIU_REAL, &pcoord));
3339566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Ncp*cdim, MPIU_REAL, &refcoord));
3349566063dSJacob Faibussowitsch     PetscCall(PetscFree(points));
335b80bf5b1SJoe Pusztay   }
3369566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
3379566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dm));
3389566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(plex, &locPhi));
3399566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(plex, &phi));
3409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vres,&vres));
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
3429566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(Vres, NULL, "-vel_res_view"));
343b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
344b80bf5b1SJoe Pusztay }
345b80bf5b1SJoe Pusztay 
346b80bf5b1SJoe Pusztay int main(int argc,char **argv)
347b80bf5b1SJoe Pusztay {
348b80bf5b1SJoe Pusztay   PetscInt          i, par;
349b80bf5b1SJoe Pusztay   PetscInt          locSize, p, d, dim, Np, step, *idx1, *idx2;
350b80bf5b1SJoe Pusztay   TS                ts;
351b80bf5b1SJoe Pusztay   DM                dm, sw;
352b80bf5b1SJoe Pusztay   AppCtx            user;
353b80bf5b1SJoe Pusztay   MPI_Comm          comm;
354b80bf5b1SJoe Pusztay   Vec               coorVec, kinVec, probVec, solution, position, momentum;
355b80bf5b1SJoe Pusztay   const PetscScalar *coorArr, *kinArr;
356b80bf5b1SJoe Pusztay   PetscReal         ftime   = 10., *probArr, *probVecArr;
357b80bf5b1SJoe Pusztay   IS                is1,is2;
358b80bf5b1SJoe Pusztay   PetscReal         *coor, *kin, *pos, *mom;
359b80bf5b1SJoe Pusztay 
3609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
361b80bf5b1SJoe Pusztay   comm = PETSC_COMM_WORLD;
3629566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
363b80bf5b1SJoe Pusztay   /* Create dm and particles */
3649566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
3659566063dSJacob Faibussowitsch   PetscCall(CreateFEM(dm, &user));
3669566063dSJacob Faibussowitsch   PetscCall(CreateParticles(dm, &sw, &user));
3679566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &user.snes));
3689566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(user.snes, dm));
3699566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
3709566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(user.snes));
371b80bf5b1SJoe Pusztay   {
372b80bf5b1SJoe Pusztay     Mat          J;
373b80bf5b1SJoe Pusztay     MatNullSpace nullSpace;
374b80bf5b1SJoe Pusztay 
3759566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &J));
3769566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace));
3779566063dSJacob Faibussowitsch     PetscCall(MatSetNullSpace(J, nullSpace));
3789566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
3799566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(user.snes, J, J, NULL, NULL));
3809566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
381b80bf5b1SJoe Pusztay   }
382b80bf5b1SJoe Pusztay   /* Place TSSolve in a loop to handle resetting the TS at every manual call of TSStep() */
3839566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
3849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
3859566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,user.stepSize));
3869566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,100000));
3879566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
388b80bf5b1SJoe Pusztay   for (step = 0; step < user.steps ; ++step){
389b80bf5b1SJoe Pusztay 
3909566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &kinVec));
3919566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &coorVec));
3929566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(kinVec, NULL, "-ic_vec_view"));
3939566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(sw, &dim));
3949566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(kinVec, &locSize));
3959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(locSize, &idx1));
3969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(locSize, &idx2));
3979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2*locSize, &probArr));
398b80bf5b1SJoe Pusztay     Np = locSize/dim;
3999566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(kinVec, &kinArr));
4009566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coorVec, &coorArr));
401b80bf5b1SJoe Pusztay     for (p=0; p<Np; ++p){
402b80bf5b1SJoe Pusztay         for (d=0; d<dim;++d) {
403b80bf5b1SJoe Pusztay             probArr[p*2*dim + d] = coorArr[p*dim+d];
404b80bf5b1SJoe Pusztay             probArr[(p*2+1)*dim + d] = kinArr[p*dim+d];
405b80bf5b1SJoe Pusztay         }
406b80bf5b1SJoe Pusztay     }
4079566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(kinVec, &kinArr));
4089566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coorVec, &coorArr));
409b80bf5b1SJoe Pusztay     /* Allocate for IS Strides that will contain x, y and vx, vy */
410b80bf5b1SJoe Pusztay     for (p = 0; p < Np; ++p) {
411b80bf5b1SJoe Pusztay       for (d = 0; d < dim; ++d) {
412b80bf5b1SJoe Pusztay         idx1[p*dim+d] = (p*2+0)*dim + d;
413b80bf5b1SJoe Pusztay         idx2[p*dim+d] = (p*2+1)*dim + d;
414b80bf5b1SJoe Pusztay       }
415b80bf5b1SJoe Pusztay     }
416b80bf5b1SJoe Pusztay 
4179566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, locSize, idx1, PETSC_OWN_POINTER, &is1));
4189566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, locSize, idx2, PETSC_OWN_POINTER, &is2));
419b80bf5b1SJoe Pusztay     /* DM needs to be set before splits so it propogates to sub TSs */
4209566063dSJacob Faibussowitsch     PetscCall(TSSetDM(ts, sw));
4219566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts,TSBASICSYMPLECTIC));
4229566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetIS(ts,"position",is1));
4239566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetIS(ts,"momentum",is2));
4249566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user));
4259566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user));
4269566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, step*user.stepSize));
427b80bf5b1SJoe Pusztay     if (step == 0) {
4289566063dSJacob Faibussowitsch       PetscCall(TSSetFromOptions(ts));
429b80bf5b1SJoe Pusztay     }
430b80bf5b1SJoe Pusztay     /* Compose vector from array for TS solve with all kinematic variables */
4319566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm,&probVec));
4329566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(probVec,1));
4339566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(probVec,PETSC_DECIDE,2*locSize));
4349566063dSJacob Faibussowitsch     PetscCall(VecSetUp(probVec));
4359566063dSJacob Faibussowitsch     PetscCall(VecGetArray(probVec,&probVecArr));
4365f80ce2aSJacob Faibussowitsch     for (i=0; i < 2*locSize; ++i) probVecArr[i] = probArr[i];
4379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(probVec,&probVecArr));
4389566063dSJacob Faibussowitsch     PetscCall(TSSetSolution(ts, probVec));
4399566063dSJacob Faibussowitsch     PetscCall(PetscFree(probArr));
4409566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(kinVec, NULL, "-ic_view"));
4419566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &kinVec));
4429566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &coorVec));
4439566063dSJacob Faibussowitsch     PetscCall(TSMonitor(ts, step, ts->ptime, ts->vec_sol));
444b80bf5b1SJoe Pusztay     if (!ts->steprollback) {
4459566063dSJacob Faibussowitsch       PetscCall(TSPreStep(ts));
446b80bf5b1SJoe Pusztay     }
4479566063dSJacob Faibussowitsch     PetscCall(TSStep(ts));
448b80bf5b1SJoe Pusztay     if (ts->steprollback) {
4499566063dSJacob Faibussowitsch       PetscCall(TSPostEvaluate(ts));
450b80bf5b1SJoe Pusztay     }
451b80bf5b1SJoe Pusztay     if (!ts->steprollback) {
452b80bf5b1SJoe Pusztay 
453b80bf5b1SJoe Pusztay       TSPostStep(ts);
4549566063dSJacob Faibussowitsch       PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coor));
4559566063dSJacob Faibussowitsch       PetscCall(DMSwarmGetField(sw, "kinematics", NULL, NULL, (void **) &kin));
4569566063dSJacob Faibussowitsch       PetscCall(TSGetSolution(ts, &solution));
4579566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(solution,is1,&position));
4589566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(solution,is2,&momentum));
4599566063dSJacob Faibussowitsch       PetscCall(VecGetArray(position, &pos));
4609566063dSJacob Faibussowitsch       PetscCall(VecGetArray(momentum, &mom));
461b80bf5b1SJoe Pusztay       for (par = 0; par < Np; ++par){
462b80bf5b1SJoe Pusztay         for (d=0; d<dim; ++d){
463b80bf5b1SJoe Pusztay           if (pos[par*dim+d] < 0.) {
464b80bf5b1SJoe Pusztay             coor[par*dim+d] = pos[par*dim+d] + 2.*PETSC_PI;
465b80bf5b1SJoe Pusztay           }
466b80bf5b1SJoe Pusztay           else if (pos[par*dim+d] > 2.*PETSC_PI) {
467b80bf5b1SJoe Pusztay             coor[par*dim+d] = pos[par*dim+d] - 2.*PETSC_PI;
468b80bf5b1SJoe Pusztay           }
469b80bf5b1SJoe Pusztay           else{
470b80bf5b1SJoe Pusztay             coor[par*dim+d] = pos[par*dim+d];
471b80bf5b1SJoe Pusztay           }
472b80bf5b1SJoe Pusztay           kin[par*dim+d] = mom[par*dim+d];
473b80bf5b1SJoe Pusztay         }
474b80bf5b1SJoe Pusztay       }
4759566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(position, &pos));
4769566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(momentum, &mom));
4779566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(solution,is1,&position));
4789566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(solution,is2,&momentum));
4799566063dSJacob Faibussowitsch       PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coor));
4809566063dSJacob Faibussowitsch       PetscCall(DMSwarmRestoreField(sw, "kinematics", NULL, NULL, (void **) &kin));
481b80bf5b1SJoe Pusztay     }
4829566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
4839566063dSJacob Faibussowitsch     PetscCall(DMLocalizeCoordinates(sw));
4849566063dSJacob Faibussowitsch     PetscCall(TSReset(ts));
4859566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&probVec));
4869566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1));
4879566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2));
488b80bf5b1SJoe Pusztay   }
4899566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&user.snes));
4909566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
4919566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
4929566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4939566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
494b122ec5aSJacob Faibussowitsch   return 0;
495b80bf5b1SJoe Pusztay }
496b80bf5b1SJoe Pusztay 
497b80bf5b1SJoe Pusztay /*TEST
498b80bf5b1SJoe Pusztay 
499b80bf5b1SJoe Pusztay    build:
500b80bf5b1SJoe Pusztay      requires: triangle !single !complex
501b80bf5b1SJoe Pusztay    test:
502b80bf5b1SJoe Pusztay      suffix: bsi1q3
503b80bf5b1SJoe Pusztay      args: -particlesPerCell 200\
504b80bf5b1SJoe Pusztay       -petscspace_degree 2\
505b80bf5b1SJoe Pusztay       -petscfe_default_quadrature_order 3\
506b80bf5b1SJoe Pusztay       -ts_basicsymplectic_type 1\
507b80bf5b1SJoe Pusztay       -pc_type svd\
508b80bf5b1SJoe Pusztay       -uniform\
509b80bf5b1SJoe Pusztay       -sigma 1.0e-8\
510b80bf5b1SJoe Pusztay       -timeScale 2.0e-14\
511b80bf5b1SJoe Pusztay       -stepSize 1.0e-2\
512b80bf5b1SJoe Pusztay       -ts_monitor_sp_swarm\
513b80bf5b1SJoe Pusztay       -steps 10\
514b80bf5b1SJoe Pusztay       -dm_view\
515b80bf5b1SJoe Pusztay       -dm_plex_simplex 0 -dm_plex_dim 2\
516b80bf5b1SJoe Pusztay       -dm_plex_box_lower 0,-1\
517b80bf5b1SJoe Pusztay       -dm_plex_box_upper 6.283185307179586,1\
518b80bf5b1SJoe Pusztay       -dm_plex_box_bd periodic,none\
519b80bf5b1SJoe Pusztay       -dm_plex_box_faces 4,1
520b80bf5b1SJoe Pusztay    test:
521b80bf5b1SJoe Pusztay      suffix: bsi2q3
522b80bf5b1SJoe Pusztay      args: -particlesPerCell 200\
523b80bf5b1SJoe Pusztay       -petscspace_degree 2\
524b80bf5b1SJoe Pusztay       -petscfe_default_quadrature_order 3\
525b80bf5b1SJoe Pusztay       -ts_basicsymplectic_type 2\
526b80bf5b1SJoe Pusztay       -pc_type svd\
527b80bf5b1SJoe Pusztay       -uniform\
528b80bf5b1SJoe Pusztay       -sigma 1.0e-8\
529b80bf5b1SJoe Pusztay       -timeScale 2.0e-14\
530b80bf5b1SJoe Pusztay       -stepSize 1.0e-2\
531b80bf5b1SJoe Pusztay       -ts_monitor_sp_swarm\
532b80bf5b1SJoe Pusztay       -steps 10\
533b80bf5b1SJoe Pusztay       -dm_view\
534b80bf5b1SJoe Pusztay       -dm_plex_simplex 0 -dm_plex_dim 2\
535b80bf5b1SJoe Pusztay       -dm_plex_box_lower 0,-1\
536b80bf5b1SJoe Pusztay       -dm_plex_box_upper 6.283185307179586,1\
537b80bf5b1SJoe Pusztay       -dm_plex_box_bd periodic,none\
538b80bf5b1SJoe Pusztay       -dm_plex_box_faces 4,1
539b80bf5b1SJoe Pusztay TEST*/
540