xref: /petsc/src/dm/impls/swarm/tests/ex5.c (revision d7462660db4346aff638eb6de48b31d27d7bb799)
1c4762a1bSJed Brown static char help[] = "Vlasov example of particles orbiting around a central massive point.\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 {
1030602db0SMatthew G. Knepley   PetscInt  dim;
11c4762a1bSJed Brown   PetscInt  particlesPerCell; /* The number of partices per cell */
12c4762a1bSJed Brown   PetscReal momentTol;        /* Tolerance for checking moment conservation */
13c4762a1bSJed Brown   PetscBool monitor;          /* Flag for using the TS monitor */
14c4762a1bSJed Brown   PetscBool error;            /* Flag for printing the error */
15c4762a1bSJed Brown   PetscInt  ostep;            /* print the energy at each ostep time steps */
16c4762a1bSJed Brown } AppCtx;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
19c4762a1bSJed Brown {
20c4762a1bSJed Brown   PetscErrorCode ierr;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   PetscFunctionBeginUser;
23c4762a1bSJed Brown   options->monitor          = PETSC_FALSE;
24c4762a1bSJed Brown   options->error            = PETSC_FALSE;
25c4762a1bSJed Brown   options->particlesPerCell = 1;
26c4762a1bSJed Brown   options->momentTol        = 100.0*PETSC_MACHINE_EPSILON;
27c4762a1bSJed Brown   options->ostep            = 100;
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Vlasov Options", "DMPLEX");CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   PetscFunctionReturn(0);
37c4762a1bSJed Brown }
38c4762a1bSJed Brown 
39c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   PetscErrorCode ierr;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   PetscFunctionBeginUser;
4430602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
4530602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
4830602db0SMatthew G. Knepley   ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr);
49c4762a1bSJed Brown   PetscFunctionReturn(0);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   DM             dm;
55c4762a1bSJed Brown   AppCtx        *user;
56c4762a1bSJed Brown   PetscRandom    rnd;
57c4762a1bSJed Brown   PetscBool      simplex;
58c4762a1bSJed Brown   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
59c4762a1bSJed Brown   PetscInt       dim, d, cStart, cEnd, c, Np, p;
60c4762a1bSJed Brown   PetscErrorCode ierr;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBeginUser;
63c4762a1bSJed Brown   ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
66c4762a1bSJed Brown 
673ec1f749SStefano Zampini   ierr = DMGetApplicationContext(dmSw, &user);CHKERRQ(ierr);
68c4762a1bSJed Brown   Np   = user->particlesPerCell;
69c4762a1bSJed Brown   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
70c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
7130602db0SMatthew G. Knepley   ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr);
74c4762a1bSJed Brown   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
75c4762a1bSJed Brown   ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
76c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
77c4762a1bSJed Brown     if (Np == 1) {
78c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
79c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
80c4762a1bSJed Brown     } else {
81c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */
82c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
83c4762a1bSJed Brown         const PetscInt n   = c*Np + p;
84c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
85c4762a1bSJed Brown 
86c4762a1bSJed Brown         for (d = 0; d < dim; ++d) {
87c4762a1bSJed Brown           ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr);
88c4762a1bSJed Brown           sum += refcoords[d];
89c4762a1bSJed Brown         }
90c4762a1bSJed Brown         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
91c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
92c4762a1bSJed Brown       }
93c4762a1bSJed Brown     }
94c4762a1bSJed Brown   }
95c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   DM             dm;
104c4762a1bSJed Brown   AppCtx        *user;
105c4762a1bSJed Brown   PetscScalar   *initialConditions;
106c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np, p;
107c4762a1bSJed Brown   PetscErrorCode ierr;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   PetscFunctionBeginUser;
1103ec1f749SStefano Zampini   ierr = DMGetApplicationContext(dmSw, &user);CHKERRQ(ierr);
111c4762a1bSJed Brown   Np   = user->particlesPerCell;
112c4762a1bSJed Brown   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr);
116c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
117c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
118c4762a1bSJed Brown       const PetscInt n = c*Np + p;
119c4762a1bSJed Brown 
120c4762a1bSJed Brown       initialConditions[(n*2 + 0)*dim + 0] = n+1;
121c4762a1bSJed Brown       initialConditions[(n*2 + 0)*dim + 1] = 0;
122c4762a1bSJed Brown       initialConditions[(n*2 + 1)*dim + 0] = 0;
123c4762a1bSJed Brown       initialConditions[(n*2 + 1)*dim + 1] = PetscSqrtReal(1000./(n+1.));
124c4762a1bSJed Brown     }
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown   ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr);
127c4762a1bSJed Brown   PetscFunctionReturn(0);
128c4762a1bSJed Brown }
129c4762a1bSJed Brown 
130c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
131c4762a1bSJed Brown {
132c4762a1bSJed Brown   PetscInt      *cellid;
133c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
134c4762a1bSJed Brown   PetscErrorCode ierr;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   PetscFunctionBeginUser;
137c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2*dim, PETSC_REAL);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
150c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
151c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
152c4762a1bSJed Brown       const PetscInt n = c*Np + p;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown       cellid[n] = c;
155c4762a1bSJed Brown     }
156c4762a1bSJed Brown   }
157c4762a1bSJed Brown   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
160c4762a1bSJed Brown   PetscFunctionReturn(0);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown /* Create particle RHS Functions for gravity with G = 1 for simplification */
164c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
165c4762a1bSJed Brown {
166c4762a1bSJed Brown   const PetscScalar *v;
167c4762a1bSJed Brown   PetscScalar       *xres;
168c4762a1bSJed Brown   PetscInt          Np, p, dim, d;
169c4762a1bSJed Brown   PetscErrorCode    ierr;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   PetscFunctionBeginUser;
17230602db0SMatthew G. Knepley   /* The DM is not currently pushed down to the splits */
17330602db0SMatthew G. Knepley   dim  = ((AppCtx *) ctx)->dim;
174c4762a1bSJed Brown   ierr = VecGetLocalSize(Xres, &Np);CHKERRQ(ierr);
175c4762a1bSJed Brown   Np  /= dim;
176c4762a1bSJed Brown   ierr = VecGetArray(Xres, &xres);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr = VecGetArrayRead(V, &v);CHKERRQ(ierr);
178c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
179c4762a1bSJed Brown      for (d = 0; d < dim; ++d) {
180c4762a1bSJed Brown        xres[p*dim+d] = v[p*dim+d];
181c4762a1bSJed Brown      }
182c4762a1bSJed Brown   }
183c4762a1bSJed Brown   ierr = VecRestoreArrayRead(V,& v);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = VecRestoreArray(Xres, &xres);CHKERRQ(ierr);
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
189c4762a1bSJed Brown {
190c4762a1bSJed Brown   const PetscScalar *x;
191c4762a1bSJed Brown   PetscScalar       *vres;
192c4762a1bSJed Brown   PetscInt          Np, p, dim, d;
193c4762a1bSJed Brown   PetscErrorCode    ierr;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   PetscFunctionBeginUser;
19630602db0SMatthew G. Knepley   /* The DM is not currently pushed down to the splits */
19730602db0SMatthew G. Knepley   dim  = ((AppCtx *) ctx)->dim;
198c4762a1bSJed Brown   ierr = VecGetLocalSize(Vres, &Np);CHKERRQ(ierr);
199c4762a1bSJed Brown   Np  /= dim;
200c4762a1bSJed Brown   ierr = VecGetArray(Vres, &vres);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
202c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
203c4762a1bSJed Brown     const PetscScalar rsqr = DMPlex_NormD_Internal(dim, &x[p*dim]);
204c4762a1bSJed Brown 
205c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
206c4762a1bSJed Brown       vres[p*dim+d] = -(1000./(p+1.)) * x[p*dim+d]/PetscSqr(rsqr);
207c4762a1bSJed Brown     }
208c4762a1bSJed Brown   }
209c4762a1bSJed Brown   ierr = VecRestoreArray(Vres, &vres);CHKERRQ(ierr);
210c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t , Vec U, Vec R, void *ctx)
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   DM                dm;
217c4762a1bSJed Brown   const PetscScalar *u;
218c4762a1bSJed Brown   PetscScalar       *r;
219c4762a1bSJed Brown   PetscInt          Np, p, dim, d;
220c4762a1bSJed Brown   PetscErrorCode    ierr;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   PetscFunctionBeginUser;
223c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
226c4762a1bSJed Brown   Np  /= 2*dim;
227c4762a1bSJed Brown   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = VecGetArray(R, &r);CHKERRQ(ierr);
229c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
230c4762a1bSJed Brown     const PetscScalar rsqr = DMPlex_NormD_Internal(dim, &u[p*2*dim]);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
233c4762a1bSJed Brown         r[p*2*dim+d]   = u[p*2*dim+d+2];
234c4762a1bSJed Brown         r[p*2*dim+d+2] = -(1000./(1.+p)) * u[p*2*dim+d]/PetscSqr(rsqr);
235c4762a1bSJed Brown     }
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = VecRestoreArray(R, &r);CHKERRQ(ierr);
239c4762a1bSJed Brown   PetscFunctionReturn(0);
240c4762a1bSJed Brown }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u)
243c4762a1bSJed Brown {
244c4762a1bSJed Brown   DM             dm;
245c4762a1bSJed Brown   AppCtx        *user;
246c4762a1bSJed Brown   PetscErrorCode ierr;
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   PetscFunctionBeginUser;
249c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
2503ec1f749SStefano Zampini   ierr = DMGetApplicationContext(dm, &user);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = SetInitialCoordinates(dm);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = SetInitialConditions(dm, u);CHKERRQ(ierr);
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
257c4762a1bSJed Brown {
258c4762a1bSJed Brown   MPI_Comm           comm;
259c4762a1bSJed Brown   DM                 sdm;
260c4762a1bSJed Brown   AppCtx            *user;
261c4762a1bSJed Brown   const PetscScalar *u, *coords;
262c4762a1bSJed Brown   PetscScalar       *e;
263c4762a1bSJed Brown   PetscReal          t;
264c4762a1bSJed Brown   PetscInt           dim, Np, p;
265c4762a1bSJed Brown   PetscErrorCode     ierr;
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   PetscFunctionBeginUser;
268c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
269c4762a1bSJed Brown   ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr);
2703ec1f749SStefano Zampini   ierr = DMGetApplicationContext(sdm, &user);CHKERRQ(ierr);
271c4762a1bSJed Brown   ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr);
272c4762a1bSJed Brown   ierr = TSGetSolveTime(ts, &t);CHKERRQ(ierr);
273c4762a1bSJed Brown   ierr = VecGetArray(E, &e);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
276c4762a1bSJed Brown   Np  /= 2*dim;
277c4762a1bSJed Brown   ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
278c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
279c4762a1bSJed Brown     const PetscScalar *x     = &u[(p*2+0)*dim];
280c4762a1bSJed Brown     const PetscScalar *v     = &u[(p*2+1)*dim];
281c4762a1bSJed Brown     const PetscReal    x0    = p+1.;
282c4762a1bSJed Brown     const PetscReal    omega = PetscSqrtReal(1000./(p+1.))/x0;
283c4762a1bSJed Brown     const PetscReal    xe[3] = { x0*PetscCosReal(omega*t),       x0*PetscSinReal(omega*t),       0.0};
284c4762a1bSJed Brown     const PetscReal    ve[3] = {-x0*omega*PetscSinReal(omega*t), x0*omega*PetscCosReal(omega*t), 0.0};
285c4762a1bSJed Brown     PetscInt           d;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
288c4762a1bSJed Brown       e[(p*2+0)*dim+d] = x[d] - xe[d];
289c4762a1bSJed Brown       e[(p*2+1)*dim+d] = v[d] - ve[d];
290c4762a1bSJed Brown     }
29130602db0SMatthew G. Knepley     if (user->error) {ierr = PetscPrintf(comm, "t %.4g: p%D error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g\n", t, p, (double) DMPlex_NormD_Internal(dim, &e[(p*2+0)*dim]), (double) DMPlex_NormD_Internal(dim, &e[(p*2+1)*dim]), (double) x[0], (double) x[1], (double) v[0], (double) v[1], (double) xe[0], (double) xe[1], (double) ve[0], (double) ve[1], 0.5*DMPlex_NormD_Internal(dim, v), (double) (0.5*(1000./(p+1))));}
292c4762a1bSJed Brown   }
293c4762a1bSJed Brown   ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
294c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
295c4762a1bSJed Brown   ierr = VecRestoreArray(E, &e);CHKERRQ(ierr);
296c4762a1bSJed Brown   PetscFunctionReturn(0);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown int main(int argc, char **argv)
300c4762a1bSJed Brown {
301c4762a1bSJed Brown   TS                 ts;
302c4762a1bSJed Brown   TSConvergedReason  reason;
303c4762a1bSJed Brown   DM                 dm, sw;
304c4762a1bSJed Brown   Vec                u;
305c4762a1bSJed Brown   IS                 is1, is2;
306c4762a1bSJed Brown   PetscInt          *idx1, *idx2;
307c4762a1bSJed Brown   MPI_Comm           comm;
308c4762a1bSJed Brown   AppCtx             user;
309c4762a1bSJed Brown   const PetscScalar *endVals;
310c4762a1bSJed Brown   PetscReal          ftime   = .1;
311c4762a1bSJed Brown   PetscInt           locSize, dim, d, Np, p, steps;
312c4762a1bSJed Brown   PetscErrorCode     ierr;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
315c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
316c4762a1bSJed Brown   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   ierr = TSCreate(comm, &ts);CHKERRQ(ierr);
323c4762a1bSJed Brown   ierr = TSSetType(ts, TSBASICSYMPLECTIC);CHKERRQ(ierr);
324c4762a1bSJed Brown   ierr = TSSetDM(ts, sw);CHKERRQ(ierr);
325c4762a1bSJed Brown   ierr = TSSetMaxTime(ts, ftime);CHKERRQ(ierr);
326c4762a1bSJed Brown   ierr = TSSetTimeStep(ts, 0.0001);CHKERRQ(ierr);
327c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts, 10);CHKERRQ(ierr);
328c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
329c4762a1bSJed Brown   ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr);
330c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user);CHKERRQ(ierr);
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   ierr = DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr);
333c4762a1bSJed Brown   ierr = DMGetDimension(sw, &dim);CHKERRQ(ierr);
334c4762a1bSJed Brown   ierr = VecGetLocalSize(u, &locSize);CHKERRQ(ierr);
335c4762a1bSJed Brown   Np   = locSize/(2*dim);
336c4762a1bSJed Brown   ierr = PetscMalloc1(locSize/2, &idx1);CHKERRQ(ierr);
337c4762a1bSJed Brown   ierr = PetscMalloc1(locSize/2, &idx2);CHKERRQ(ierr);
338c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
339c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
340c4762a1bSJed Brown       idx1[p*dim+d] = (p*2+0)*dim + d;
341c4762a1bSJed Brown       idx2[p*dim+d] = (p*2+1)*dim + d;
342c4762a1bSJed Brown     }
343c4762a1bSJed Brown   }
344c4762a1bSJed Brown   ierr = ISCreateGeneral(comm, locSize/2, idx1, PETSC_OWN_POINTER, &is1);CHKERRQ(ierr);
345c4762a1bSJed Brown   ierr = ISCreateGeneral(comm, locSize/2, idx2, PETSC_OWN_POINTER, &is2);CHKERRQ(ierr);
346c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts, "position", is1);CHKERRQ(ierr);
347c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts, "momentum", is2);CHKERRQ(ierr);
348c4762a1bSJed Brown   ierr = ISDestroy(&is1);CHKERRQ(ierr);
349c4762a1bSJed Brown   ierr = ISDestroy(&is2);CHKERRQ(ierr);
350c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user);CHKERRQ(ierr);
351c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user);CHKERRQ(ierr);
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
354c4762a1bSJed Brown   ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr);
355c4762a1bSJed Brown   ierr = TSSetComputeExactError(ts, ComputeError);CHKERRQ(ierr);
356c4762a1bSJed Brown   ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr);
357c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-init_view");CHKERRQ(ierr);
358c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
359c4762a1bSJed Brown   ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr);
360c4762a1bSJed Brown   ierr = TSGetConvergedReason(ts, &reason);CHKERRQ(ierr);
361c4762a1bSJed Brown   ierr = TSGetStepNumber(ts, &steps);CHKERRQ(ierr);
362c4762a1bSJed Brown   ierr = PetscPrintf(comm,"%s at time %g after %D steps\n", TSConvergedReasons[reason], (double) ftime, steps);CHKERRQ(ierr);
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   ierr = VecGetArrayRead(u, &endVals);CHKERRQ(ierr);
365c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
366c4762a1bSJed Brown     const PetscReal norm = DMPlex_NormD_Internal(dim, &endVals[(p*2 + 1)*dim]);
367c4762a1bSJed Brown     ierr = PetscPrintf(comm, "Particle %D initial Energy: %g  Final Energy: %g\n", p, (double) (0.5*(1000./(p+1))), (double) 0.5*PetscSqr(norm));CHKERRQ(ierr);
368c4762a1bSJed Brown   }
369c4762a1bSJed Brown   ierr = VecRestoreArrayRead(u, &endVals);CHKERRQ(ierr);
370c4762a1bSJed Brown   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr);
371c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
372c4762a1bSJed Brown   ierr = DMDestroy(&sw);CHKERRQ(ierr);
373c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
374c4762a1bSJed Brown   ierr = PetscFinalize();
375c4762a1bSJed Brown   return ierr;
376c4762a1bSJed Brown }
377c4762a1bSJed Brown 
378c4762a1bSJed Brown /*TEST
379c4762a1bSJed Brown 
380c4762a1bSJed Brown    build:
381c4762a1bSJed Brown      requires: triangle !single !complex
382c4762a1bSJed Brown    test:
383c4762a1bSJed Brown      suffix: bsi1
384*d7462660SMatthew Knepley      args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
385*d7462660SMatthew Knepley            -ts_basicsymplectic_type 1 -ts_max_time 0.1 -ts_convergence_estimate -convest_num_refine 2 \
386*d7462660SMatthew Knepley            -dm_view -sw_view -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -ts_monitor_sp_swarm_phase 0
387c4762a1bSJed Brown    test:
388c4762a1bSJed Brown      suffix: bsi2
389*d7462660SMatthew Knepley      args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
390*d7462660SMatthew Knepley            -ts_basicsymplectic_type 2 -ts_max_time 0.1 -ts_convergence_estimate -convest_num_refine 2 \
391*d7462660SMatthew Knepley            -dm_view -sw_view -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -ts_monitor_sp_swarm_phase 0
392c4762a1bSJed Brown    test:
393c4762a1bSJed Brown      suffix: euler
394*d7462660SMatthew Knepley      args: -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
395*d7462660SMatthew Knepley            -ts_type euler -ts_max_time 0.1 -ts_convergence_estimate -convest_num_refine 2 \
396*d7462660SMatthew Knepley            -dm_view -sw_view -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -ts_monitor_sp_swarm_phase 0
397c4762a1bSJed Brown 
398c4762a1bSJed Brown TEST*/
399