xref: /petsc/src/ts/tutorials/hamiltonian/ex2.c (revision 9566063d113dddea24716c546802770db7481bc0)
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   PetscErrorCode ierr;
35b80bf5b1SJoe Pusztay 
36b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
37b80bf5b1SJoe Pusztay   options->dim              = 2;
38b80bf5b1SJoe Pusztay   options->simplex          = PETSC_TRUE;
39b80bf5b1SJoe Pusztay   options->monitor          = PETSC_TRUE;
40b80bf5b1SJoe Pusztay   options->particlesPerCell = 1;
41b80bf5b1SJoe Pusztay   options->k                = 1;
42b80bf5b1SJoe Pusztay   options->particleRelDx    = 1.e-20;
43b80bf5b1SJoe Pusztay   options->momentTol        = 100.*PETSC_MACHINE_EPSILON;
44b80bf5b1SJoe Pusztay   options->sigma            = 1.;
45b80bf5b1SJoe Pusztay   options->timeScale        = 1.0e-6;
46b80bf5b1SJoe Pusztay   options->uniform          = PETSC_FALSE;
47b80bf5b1SJoe Pusztay   options->steps            = 1;
48b80bf5b1SJoe Pusztay   options->stepSize         = 0.01;
49b80bf5b1SJoe Pusztay   options->bdm              = PETSC_FALSE;
50b80bf5b1SJoe Pusztay 
51*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Two Stream options", "DMPLEX");PetscCall(ierr);
52*9566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(options->meshFilename, ""));
53*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-dim", "The topological mesh dimension", "ex2.c", options->dim, &options->dim, NULL));
54*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-steps", "TS steps to take", "ex2.c", options->steps, &options->steps, NULL));
55*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-monitor", "To use the TS monitor or not", "ex2.c", options->monitor, &options->monitor, NULL));
56*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex2.c", options->simplex, &options->simplex, NULL));
57*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-uniform", "Uniform particle spacing", "ex2.c", options->uniform, &options->uniform, NULL));
58*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-bdm", "Use H1 instead of C0", "ex2.c", options->bdm, &options->bdm, NULL));
59*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex2.c", options->meshFilename, options->meshFilename, PETSC_MAX_PATH_LEN, NULL));
60*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex5.c", options->k, &options->k, NULL));
61*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL));
62*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-sigma","parameter","<1>",options->sigma,&options->sigma,PETSC_NULL));
63*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-stepSize","parameter","<1e-2>",options->stepSize,&options->stepSize,PETSC_NULL));
64*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-timeScale","parameter","<1>",options->timeScale,&options->timeScale,PETSC_NULL));
65*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL));
66*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
67b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
68b80bf5b1SJoe Pusztay }
69b80bf5b1SJoe Pusztay 
70b80bf5b1SJoe Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
71b80bf5b1SJoe Pusztay {
72b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
73*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
74*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
75*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
76*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
77b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
78b80bf5b1SJoe Pusztay }
79b80bf5b1SJoe Pusztay 
80b80bf5b1SJoe Pusztay static void laplacian_f1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
81b80bf5b1SJoe Pusztay                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
82b80bf5b1SJoe Pusztay                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
83b80bf5b1SJoe Pusztay                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
84b80bf5b1SJoe Pusztay {
85b80bf5b1SJoe Pusztay   PetscInt d;
86b80bf5b1SJoe Pusztay   for (d = 0; d < dim; ++d) {f1[d] = u_x[d];}
87b80bf5b1SJoe Pusztay }
88b80bf5b1SJoe Pusztay 
89b80bf5b1SJoe Pusztay static void laplacian(PetscInt dim, PetscInt Nf, PetscInt NfAux,
90b80bf5b1SJoe Pusztay                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
91b80bf5b1SJoe Pusztay                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
92b80bf5b1SJoe Pusztay                       PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
93b80bf5b1SJoe Pusztay {
94b80bf5b1SJoe Pusztay   PetscInt d;
95b80bf5b1SJoe Pusztay   for (d = 0; d < dim; ++d) {g3[d*dim+d] = 1.0;}
96b80bf5b1SJoe Pusztay }
97b80bf5b1SJoe Pusztay 
98b80bf5b1SJoe Pusztay static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
99b80bf5b1SJoe Pusztay {
100b80bf5b1SJoe Pusztay   PetscFE        fe;
101b80bf5b1SJoe Pusztay   PetscDS        ds;
102b80bf5b1SJoe Pusztay   DMPolytopeType ct;
103b80bf5b1SJoe Pusztay   PetscBool      simplex;
104b80bf5b1SJoe Pusztay   PetscInt       dim, cStart;
105b80bf5b1SJoe Pusztay 
106b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
107*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
108*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
109*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
110b80bf5b1SJoe Pusztay   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
111*9566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, NULL, -1, &fe));
112*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, "potential"));
113*9566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
114*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
115*9566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
116*9566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
117*9566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, NULL, laplacian_f1));
118*9566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, laplacian));
119b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
120b80bf5b1SJoe Pusztay }
121b80bf5b1SJoe Pusztay 
122b80bf5b1SJoe Pusztay /*
123b80bf5b1SJoe Pusztay   Initialize particle coordinates uniformly and with opposing velocities
124b80bf5b1SJoe Pusztay */
125b80bf5b1SJoe Pusztay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
126b80bf5b1SJoe Pusztay {
127b80bf5b1SJoe Pusztay   PetscRandom    rnd, rndp;
128b80bf5b1SJoe Pusztay   PetscReal      interval = user->particleRelDx;
129b80bf5b1SJoe Pusztay   PetscScalar    value, *vals;
130b80bf5b1SJoe Pusztay   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ, *initialConditions, normalized_vel;
131b80bf5b1SJoe Pusztay   PetscInt      *cellid, cStart;
132b80bf5b1SJoe Pusztay   PetscInt       Ncell, Np = user->particlesPerCell, p, c, dim, d;
133b80bf5b1SJoe Pusztay 
134b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
135*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
136*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
137*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
138*9566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
139*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rnd));
140*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rnd, 0.0, 1.0));
141*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rnd));
142*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject) dm), &rndp));
143*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rndp, -interval, interval));
144*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rndp));
145*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
146*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
147*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
148*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", dim, PETSC_REAL));
149*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
150*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &Ncell));
151*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0));
152*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
153*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
154*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
155*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **) &vals));
156*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(*sw, "kinematics", NULL, NULL, (void **) &initialConditions));
157*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ));
158b80bf5b1SJoe Pusztay   for (c = cStart; c < Ncell; c++) {
159b80bf5b1SJoe Pusztay     if (Np == 1) {
160*9566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
161b80bf5b1SJoe Pusztay       cellid[c] = c;
162b80bf5b1SJoe Pusztay       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
163b80bf5b1SJoe Pusztay     } else {
164b80bf5b1SJoe Pusztay       for (d = 0; d < dim; ++d) xi0[d] = -1.0;
165*9566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */
166b80bf5b1SJoe Pusztay       for (p = 0; p < Np; ++p) {
167b80bf5b1SJoe Pusztay         const PetscInt n   = c*Np + p;
168b80bf5b1SJoe Pusztay         PetscReal      refcoords[3], spacing;
169b80bf5b1SJoe Pusztay 
170b80bf5b1SJoe Pusztay         cellid[n] = c;
171b80bf5b1SJoe Pusztay         if (user->uniform) {
172b80bf5b1SJoe Pusztay           spacing = 2./Np;
173*9566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValue(rnd, &value));
174b80bf5b1SJoe Pusztay           for (d=0; d<dim; ++d) refcoords[d] = d == 0 ? -1. + spacing/2. + p*spacing + value/100. : 0.;
175b80bf5b1SJoe Pusztay         }
176b80bf5b1SJoe Pusztay         else{
177*9566063dSJacob Faibussowitsch           for (d = 0; d < dim; ++d) {PetscCall(PetscRandomGetValue(rnd, &value)); refcoords[d] = d == 0 ? PetscRealPart(value) : 0. ;}
178b80bf5b1SJoe Pusztay         }
179b80bf5b1SJoe Pusztay         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
180b80bf5b1SJoe Pusztay         /* constant particle weights */
181b80bf5b1SJoe Pusztay         for (d = 0; d < dim; ++d) vals[n] = user->sigma/Np;
182b80bf5b1SJoe Pusztay       }
183b80bf5b1SJoe Pusztay     }
184b80bf5b1SJoe Pusztay   }
185*9566063dSJacob Faibussowitsch   PetscCall(PetscFree5(centroid, xi0, v0, J, invJ));
186b80bf5b1SJoe Pusztay   normalized_vel = 1.;
187b80bf5b1SJoe Pusztay   for (c = 0; c < Ncell; ++c) {
188b80bf5b1SJoe Pusztay     for (p = 0; p < Np; ++p) {
189b80bf5b1SJoe Pusztay       if (p%2 == 0) {
190b80bf5b1SJoe Pusztay         for (d = 0; d < dim; ++d) initialConditions[(c*Np + p)*dim + d] = d == 0 ? normalized_vel : 0.;
191b80bf5b1SJoe Pusztay       }
192b80bf5b1SJoe Pusztay       else {
193b80bf5b1SJoe Pusztay         for (d = 0; d < dim; ++d) initialConditions[(c*Np + p)*dim + d] = d == 0 ? -(normalized_vel) : 0.;
194b80bf5b1SJoe Pusztay       }
195b80bf5b1SJoe Pusztay     }
196b80bf5b1SJoe Pusztay   }
197*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
198*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid));
199*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **) &vals));
200*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(*sw, "kinematics", NULL, NULL, (void **) &initialConditions));
201*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rnd));
202*9566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rndp));
203*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
204*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
205*9566063dSJacob Faibussowitsch   PetscCall(DMLocalizeCoordinates(*sw));
206b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
207b80bf5b1SJoe Pusztay }
208b80bf5b1SJoe Pusztay 
209b80bf5b1SJoe Pusztay /* Solve for particle position updates */
210b80bf5b1SJoe Pusztay static PetscErrorCode RHSFunction1(TS ts,PetscReal t,Vec V,Vec Posres,void *ctx)
211b80bf5b1SJoe Pusztay {
212b80bf5b1SJoe Pusztay   const PetscScalar *v;
213b80bf5b1SJoe Pusztay   PetscScalar       *posres;
214b80bf5b1SJoe Pusztay   PetscInt          Np, p, dim, d;
215b80bf5b1SJoe Pusztay   DM                dm;
216b80bf5b1SJoe Pusztay 
217b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
218*9566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Posres, &Np));
219*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Posres,&posres));
220*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V,&v));
221*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
222*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
223b80bf5b1SJoe Pusztay   Np  /= dim;
224b80bf5b1SJoe Pusztay   for (p = 0; p < Np; ++p) {
225b80bf5b1SJoe Pusztay      for (d = 0; d < dim; ++d) {
226b80bf5b1SJoe Pusztay        posres[p*dim+d] = v[p*dim+d];
227b80bf5b1SJoe Pusztay      }
228b80bf5b1SJoe Pusztay   }
229*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V,&v));
230*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Posres,&posres));
231b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
232b80bf5b1SJoe Pusztay 
233b80bf5b1SJoe Pusztay }
234b80bf5b1SJoe Pusztay 
235b80bf5b1SJoe Pusztay /*
236b80bf5b1SJoe Pusztay   Solve for the gradient of the electric field and apply force to particles.
237b80bf5b1SJoe Pusztay  */
238b80bf5b1SJoe Pusztay static PetscErrorCode RHSFunction2(TS ts,PetscReal t,Vec X,Vec Vres,void *ctx)
239b80bf5b1SJoe Pusztay {
240b80bf5b1SJoe Pusztay  AppCtx             *user = (AppCtx *) ctx;
241b80bf5b1SJoe Pusztay   DM                 dm, plex;
242b80bf5b1SJoe Pusztay   PetscDS            prob;
243b80bf5b1SJoe Pusztay   PetscFE            fe;
244b80bf5b1SJoe Pusztay   Mat                M_p;
245b80bf5b1SJoe Pusztay   Vec                phi, locPhi, rho, f;
246b80bf5b1SJoe Pusztay   const PetscScalar *x;
247b80bf5b1SJoe Pusztay   PetscScalar       *vres;
248b80bf5b1SJoe Pusztay   PetscReal         *coords, phi_0;
249b80bf5b1SJoe Pusztay   PetscInt           dim, d, cStart, cEnd, cell, cdim;
250b80bf5b1SJoe Pusztay   PetscReal          m_e = 9.11e-31, q_e = 1.60e-19, epsi_0 = 8.85e-12;
251b80bf5b1SJoe Pusztay 
252b80bf5b1SJoe Pusztay   PetscFunctionBeginUser;
253b80bf5b1SJoe Pusztay   PetscObjectSetName((PetscObject) X, "rhsf2 position");
254b80bf5b1SJoe Pusztay   VecViewFromOptions(X, NULL, "-rhsf2_x_view");
255*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
256*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vres,&vres));
257*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
258*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
259*9566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(user->snes, &plex));
260*9566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(plex, &cdim));
261*9566063dSJacob Faibussowitsch   PetscCall(DMGetDS(plex, &prob));
262*9566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe));
263*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(plex, &phi));
264*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(plex, &locPhi));
265*9566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dm, plex, &M_p));
266*9566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
267*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(plex, &rho));
268*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dm, "w_q", &f));
269*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) f, "weights vector"));
270*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
271*9566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(M_p, f, rho));
272*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, "w_q", &f));
273*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) rho, "rho"));
274*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
275b80bf5b1SJoe Pusztay   /* Take nullspace out of rhs */
276b80bf5b1SJoe Pusztay   {
277b80bf5b1SJoe Pusztay     PetscScalar sum;
278b80bf5b1SJoe Pusztay     PetscInt    n;
279b80bf5b1SJoe Pusztay     phi_0 = (user->sigma*user->sigma*user->sigma)*(user->timeScale*user->timeScale)/(m_e*q_e*epsi_0);
280b80bf5b1SJoe Pusztay 
281*9566063dSJacob Faibussowitsch     PetscCall(VecGetSize(rho, &n));
282*9566063dSJacob Faibussowitsch     PetscCall(VecSum(rho, &sum));
283*9566063dSJacob Faibussowitsch     PetscCall(VecShift(rho, -sum/n));
284b80bf5b1SJoe Pusztay 
285*9566063dSJacob Faibussowitsch     PetscCall(VecSum(rho, &sum));
2863c633725SBarry Smith     PetscCheck(PetscAbsScalar(sum) <= 1.0e-10,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Charge should have no DC component %g", sum);
287*9566063dSJacob Faibussowitsch     PetscCall(VecScale(rho, phi_0));
288b80bf5b1SJoe Pusztay   }
289*9566063dSJacob Faibussowitsch   PetscCall(VecSet(phi, 0.0));
290*9566063dSJacob Faibussowitsch   PetscCall(SNESSolve(user->snes, rho, phi));
291*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
292*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(plex, &rho));
293*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&M_p));
294*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(plex, phi, INSERT_VALUES, locPhi));
295*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(plex, phi, INSERT_VALUES, locPhi));
296*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortGetAccess(dm));
297*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
298*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
299b80bf5b1SJoe Pusztay   for (cell = cStart; cell < cEnd; ++cell) {
300b80bf5b1SJoe Pusztay     PetscTabulation tab;
301b80bf5b1SJoe Pusztay     PetscReal    v[3], J[9], invJ[9], detJ;
302b80bf5b1SJoe Pusztay     PetscScalar *ph       = PETSC_NULL;
303b80bf5b1SJoe Pusztay     PetscReal   *pcoord   = PETSC_NULL;
304b80bf5b1SJoe Pusztay     PetscReal   *refcoord = PETSC_NULL;
305b80bf5b1SJoe Pusztay     PetscInt    *points   = PETSC_NULL, Ncp, cp;
306b80bf5b1SJoe Pusztay     PetscScalar  gradPhi[3];
307b80bf5b1SJoe Pusztay 
308*9566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, NULL, v, J, invJ, &detJ));
309*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmSortGetPointsPerCell(dm, cell, &Ncp, &points));
310*9566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, Ncp*cdim, MPIU_REAL, &pcoord));
311*9566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(dm, Ncp*cdim, MPIU_REAL, &refcoord));
312b80bf5b1SJoe Pusztay     for (cp = 0; cp < Ncp; ++cp) {
313b80bf5b1SJoe Pusztay       for (d = 0; d < cdim; ++d) {
314b80bf5b1SJoe Pusztay         pcoord[cp*cdim+d] = coords[points[cp]*cdim+d];
315b80bf5b1SJoe Pusztay       }
316b80bf5b1SJoe Pusztay     }
317*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCoordinatesToReference(plex, cell, Ncp, pcoord, refcoord));
318*9566063dSJacob Faibussowitsch     PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
319*9566063dSJacob Faibussowitsch     PetscCall(DMPlexVecGetClosure(plex, NULL, locPhi, cell, NULL, &ph));
320b80bf5b1SJoe Pusztay     for (cp = 0; cp < Ncp; ++cp) {
321b80bf5b1SJoe Pusztay       const PetscInt p = points[cp];
322b80bf5b1SJoe Pusztay       gradPhi[0] = 0.0;
323b80bf5b1SJoe Pusztay       gradPhi[1] = 0.0;
324b80bf5b1SJoe Pusztay       gradPhi[2] = 0.0;
325b80bf5b1SJoe Pusztay       const PetscReal *basisDer = tab->T[1];
326b80bf5b1SJoe Pusztay 
327*9566063dSJacob Faibussowitsch       PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, ph, cdim, invJ, NULL, cp, gradPhi));
328b80bf5b1SJoe Pusztay       for (d = 0; d < cdim; ++d) {
329b80bf5b1SJoe Pusztay         vres[p*cdim+d] = d == 0 ? gradPhi[d] : 0.;
330b80bf5b1SJoe Pusztay       }
331b80bf5b1SJoe Pusztay     }
332*9566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(plex, NULL, locPhi, cell, NULL, &ph));
333*9566063dSJacob Faibussowitsch     PetscCall(PetscTabulationDestroy(&tab));
334*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Ncp*cdim, MPIU_REAL, &pcoord));
335*9566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(dm, Ncp*cdim, MPIU_REAL, &refcoord));
336*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(points));
337b80bf5b1SJoe Pusztay   }
338*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords));
339*9566063dSJacob Faibussowitsch   PetscCall(DMSwarmSortRestoreAccess(dm));
340*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(plex, &locPhi));
341*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(plex, &phi));
342*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vres,&vres));
343*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
344*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(Vres, NULL, "-vel_res_view"));
345b80bf5b1SJoe Pusztay   PetscFunctionReturn(0);
346b80bf5b1SJoe Pusztay }
347b80bf5b1SJoe Pusztay 
348b80bf5b1SJoe Pusztay int main(int argc,char **argv)
349b80bf5b1SJoe Pusztay {
350b80bf5b1SJoe Pusztay   PetscInt          i, par;
351b80bf5b1SJoe Pusztay   PetscInt          locSize, p, d, dim, Np, step, *idx1, *idx2;
352b80bf5b1SJoe Pusztay   TS                ts;
353b80bf5b1SJoe Pusztay   DM                dm, sw;
354b80bf5b1SJoe Pusztay   AppCtx            user;
355b80bf5b1SJoe Pusztay   MPI_Comm          comm;
356b80bf5b1SJoe Pusztay   Vec               coorVec, kinVec, probVec, solution, position, momentum;
357b80bf5b1SJoe Pusztay   const PetscScalar *coorArr, *kinArr;
358b80bf5b1SJoe Pusztay   PetscReal         ftime   = 10., *probArr, *probVecArr;
359b80bf5b1SJoe Pusztay   IS                is1,is2;
360b80bf5b1SJoe Pusztay   PetscReal         *coor, *kin, *pos, *mom;
361b80bf5b1SJoe Pusztay 
362*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
363b80bf5b1SJoe Pusztay   comm = PETSC_COMM_WORLD;
364*9566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
365b80bf5b1SJoe Pusztay   /* Create dm and particles */
366*9566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm, &user));
367*9566063dSJacob Faibussowitsch   PetscCall(CreateFEM(dm, &user));
368*9566063dSJacob Faibussowitsch   PetscCall(CreateParticles(dm, &sw, &user));
369*9566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &user.snes));
370*9566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(user.snes, dm));
371*9566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
372*9566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(user.snes));
373b80bf5b1SJoe Pusztay   {
374b80bf5b1SJoe Pusztay     Mat          J;
375b80bf5b1SJoe Pusztay     MatNullSpace nullSpace;
376b80bf5b1SJoe Pusztay 
377*9566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(dm, &J));
378*9566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_TRUE, 0, NULL, &nullSpace));
379*9566063dSJacob Faibussowitsch     PetscCall(MatSetNullSpace(J, nullSpace));
380*9566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
381*9566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(user.snes, J, J, NULL, NULL));
382*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
383b80bf5b1SJoe Pusztay   }
384b80bf5b1SJoe Pusztay   /* Place TSSolve in a loop to handle resetting the TS at every manual call of TSStep() */
385*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(comm, &ts));
386*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
387*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,user.stepSize));
388*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,100000));
389*9566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
390b80bf5b1SJoe Pusztay   for (step = 0; step < user.steps ; ++step){
391b80bf5b1SJoe Pusztay 
392*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &kinVec));
393*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &coorVec));
394*9566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(kinVec, NULL, "-ic_vec_view"));
395*9566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(sw, &dim));
396*9566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(kinVec, &locSize));
397*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(locSize, &idx1));
398*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(locSize, &idx2));
399*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(2*locSize, &probArr));
400b80bf5b1SJoe Pusztay     Np = locSize/dim;
401*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(kinVec, &kinArr));
402*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coorVec, &coorArr));
403b80bf5b1SJoe Pusztay     for (p=0; p<Np; ++p){
404b80bf5b1SJoe Pusztay         for (d=0; d<dim;++d) {
405b80bf5b1SJoe Pusztay             probArr[p*2*dim + d] = coorArr[p*dim+d];
406b80bf5b1SJoe Pusztay             probArr[(p*2+1)*dim + d] = kinArr[p*dim+d];
407b80bf5b1SJoe Pusztay         }
408b80bf5b1SJoe Pusztay     }
409*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(kinVec, &kinArr));
410*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coorVec, &coorArr));
411b80bf5b1SJoe Pusztay     /* Allocate for IS Strides that will contain x, y and vx, vy */
412b80bf5b1SJoe Pusztay     for (p = 0; p < Np; ++p) {
413b80bf5b1SJoe Pusztay       for (d = 0; d < dim; ++d) {
414b80bf5b1SJoe Pusztay         idx1[p*dim+d] = (p*2+0)*dim + d;
415b80bf5b1SJoe Pusztay         idx2[p*dim+d] = (p*2+1)*dim + d;
416b80bf5b1SJoe Pusztay       }
417b80bf5b1SJoe Pusztay     }
418b80bf5b1SJoe Pusztay 
419*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, locSize, idx1, PETSC_OWN_POINTER, &is1));
420*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, locSize, idx2, PETSC_OWN_POINTER, &is2));
421b80bf5b1SJoe Pusztay     /* DM needs to be set before splits so it propogates to sub TSs */
422*9566063dSJacob Faibussowitsch     PetscCall(TSSetDM(ts, sw));
423*9566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts,TSBASICSYMPLECTIC));
424*9566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetIS(ts,"position",is1));
425*9566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetIS(ts,"momentum",is2));
426*9566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user));
427*9566063dSJacob Faibussowitsch     PetscCall(TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user));
428*9566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, step*user.stepSize));
429b80bf5b1SJoe Pusztay     if (step == 0) {
430*9566063dSJacob Faibussowitsch       PetscCall(TSSetFromOptions(ts));
431b80bf5b1SJoe Pusztay     }
432b80bf5b1SJoe Pusztay     /* Compose vector from array for TS solve with all kinematic variables */
433*9566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm,&probVec));
434*9566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(probVec,1));
435*9566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(probVec,PETSC_DECIDE,2*locSize));
436*9566063dSJacob Faibussowitsch     PetscCall(VecSetUp(probVec));
437*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(probVec,&probVecArr));
4385f80ce2aSJacob Faibussowitsch     for (i=0; i < 2*locSize; ++i) probVecArr[i] = probArr[i];
439*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(probVec,&probVecArr));
440*9566063dSJacob Faibussowitsch     PetscCall(TSSetSolution(ts, probVec));
441*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(probArr));
442*9566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(kinVec, NULL, "-ic_view"));
443*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &kinVec));
444*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &coorVec));
445*9566063dSJacob Faibussowitsch     PetscCall(TSMonitor(ts, step, ts->ptime, ts->vec_sol));
446b80bf5b1SJoe Pusztay     if (!ts->steprollback) {
447*9566063dSJacob Faibussowitsch       PetscCall(TSPreStep(ts));
448b80bf5b1SJoe Pusztay     }
449*9566063dSJacob Faibussowitsch     PetscCall(TSStep(ts));
450b80bf5b1SJoe Pusztay     if (ts->steprollback) {
451*9566063dSJacob Faibussowitsch       PetscCall(TSPostEvaluate(ts));
452b80bf5b1SJoe Pusztay     }
453b80bf5b1SJoe Pusztay     if (!ts->steprollback) {
454b80bf5b1SJoe Pusztay 
455b80bf5b1SJoe Pusztay       TSPostStep(ts);
456*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coor));
457*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmGetField(sw, "kinematics", NULL, NULL, (void **) &kin));
458*9566063dSJacob Faibussowitsch       PetscCall(TSGetSolution(ts, &solution));
459*9566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(solution,is1,&position));
460*9566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(solution,is2,&momentum));
461*9566063dSJacob Faibussowitsch       PetscCall(VecGetArray(position, &pos));
462*9566063dSJacob Faibussowitsch       PetscCall(VecGetArray(momentum, &mom));
463b80bf5b1SJoe Pusztay       for (par = 0; par < Np; ++par){
464b80bf5b1SJoe Pusztay         for (d=0; d<dim; ++d){
465b80bf5b1SJoe Pusztay           if (pos[par*dim+d] < 0.) {
466b80bf5b1SJoe Pusztay             coor[par*dim+d] = pos[par*dim+d] + 2.*PETSC_PI;
467b80bf5b1SJoe Pusztay           }
468b80bf5b1SJoe Pusztay           else if (pos[par*dim+d] > 2.*PETSC_PI) {
469b80bf5b1SJoe Pusztay             coor[par*dim+d] = pos[par*dim+d] - 2.*PETSC_PI;
470b80bf5b1SJoe Pusztay           }
471b80bf5b1SJoe Pusztay           else{
472b80bf5b1SJoe Pusztay             coor[par*dim+d] = pos[par*dim+d];
473b80bf5b1SJoe Pusztay           }
474b80bf5b1SJoe Pusztay           kin[par*dim+d] = mom[par*dim+d];
475b80bf5b1SJoe Pusztay         }
476b80bf5b1SJoe Pusztay       }
477*9566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(position, &pos));
478*9566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(momentum, &mom));
479*9566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(solution,is1,&position));
480*9566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(solution,is2,&momentum));
481*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **) &coor));
482*9566063dSJacob Faibussowitsch       PetscCall(DMSwarmRestoreField(sw, "kinematics", NULL, NULL, (void **) &kin));
483b80bf5b1SJoe Pusztay     }
484*9566063dSJacob Faibussowitsch     PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
485*9566063dSJacob Faibussowitsch     PetscCall(DMLocalizeCoordinates(sw));
486*9566063dSJacob Faibussowitsch     PetscCall(TSReset(ts));
487*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&probVec));
488*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1));
489*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2));
490b80bf5b1SJoe Pusztay   }
491*9566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&user.snes));
492*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
493*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
494*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
495*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
496b122ec5aSJacob Faibussowitsch   return 0;
497b80bf5b1SJoe Pusztay }
498b80bf5b1SJoe Pusztay 
499b80bf5b1SJoe Pusztay /*TEST
500b80bf5b1SJoe Pusztay 
501b80bf5b1SJoe Pusztay    build:
502b80bf5b1SJoe Pusztay      requires: triangle !single !complex
503b80bf5b1SJoe Pusztay    test:
504b80bf5b1SJoe Pusztay      suffix: bsi1q3
505b80bf5b1SJoe Pusztay      args: -particlesPerCell 200\
506b80bf5b1SJoe Pusztay       -petscspace_degree 2\
507b80bf5b1SJoe Pusztay       -petscfe_default_quadrature_order 3\
508b80bf5b1SJoe Pusztay       -ts_basicsymplectic_type 1\
509b80bf5b1SJoe Pusztay       -pc_type svd\
510b80bf5b1SJoe Pusztay       -uniform\
511b80bf5b1SJoe Pusztay       -sigma 1.0e-8\
512b80bf5b1SJoe Pusztay       -timeScale 2.0e-14\
513b80bf5b1SJoe Pusztay       -stepSize 1.0e-2\
514b80bf5b1SJoe Pusztay       -ts_monitor_sp_swarm\
515b80bf5b1SJoe Pusztay       -steps 10\
516b80bf5b1SJoe Pusztay       -dm_view\
517b80bf5b1SJoe Pusztay       -dm_plex_simplex 0 -dm_plex_dim 2\
518b80bf5b1SJoe Pusztay       -dm_plex_box_lower 0,-1\
519b80bf5b1SJoe Pusztay       -dm_plex_box_upper 6.283185307179586,1\
520b80bf5b1SJoe Pusztay       -dm_plex_box_bd periodic,none\
521b80bf5b1SJoe Pusztay       -dm_plex_box_faces 4,1
522b80bf5b1SJoe Pusztay    test:
523b80bf5b1SJoe Pusztay      suffix: bsi2q3
524b80bf5b1SJoe Pusztay      args: -particlesPerCell 200\
525b80bf5b1SJoe Pusztay       -petscspace_degree 2\
526b80bf5b1SJoe Pusztay       -petscfe_default_quadrature_order 3\
527b80bf5b1SJoe Pusztay       -ts_basicsymplectic_type 2\
528b80bf5b1SJoe Pusztay       -pc_type svd\
529b80bf5b1SJoe Pusztay       -uniform\
530b80bf5b1SJoe Pusztay       -sigma 1.0e-8\
531b80bf5b1SJoe Pusztay       -timeScale 2.0e-14\
532b80bf5b1SJoe Pusztay       -stepSize 1.0e-2\
533b80bf5b1SJoe Pusztay       -ts_monitor_sp_swarm\
534b80bf5b1SJoe Pusztay       -steps 10\
535b80bf5b1SJoe Pusztay       -dm_view\
536b80bf5b1SJoe Pusztay       -dm_plex_simplex 0 -dm_plex_dim 2\
537b80bf5b1SJoe Pusztay       -dm_plex_box_lower 0,-1\
538b80bf5b1SJoe Pusztay       -dm_plex_box_upper 6.283185307179586,1\
539b80bf5b1SJoe Pusztay       -dm_plex_box_bd periodic,none\
540b80bf5b1SJoe Pusztay       -dm_plex_box_faces 4,1
541b80bf5b1SJoe Pusztay TEST*/
542