xref: /petsc/src/dm/impls/swarm/tests/ex4.c (revision 362febeeeb69b91ebadcb4b2dc0a22cb6dfc4097)
1c4762a1bSJed Brown static char help[] = "Example of simple hamiltonian system (harmonic oscillator) with particles and a basic symplectic integrator\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   PetscReal omega;            /* Oscillation frequency omega */
12c4762a1bSJed Brown   PetscInt  particlesPerCell; /* The number of partices per cell */
13c4762a1bSJed Brown   PetscReal momentTol;        /* Tolerance for checking moment conservation */
14c4762a1bSJed Brown   PetscBool monitor;          /* Flag for using the TS monitor */
15c4762a1bSJed Brown   PetscBool error;            /* Flag for printing the error */
16c4762a1bSJed Brown   PetscInt  ostep;            /* print the energy at each ostep time steps */
17c4762a1bSJed Brown } AppCtx;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
20c4762a1bSJed Brown {
21c4762a1bSJed Brown   PetscErrorCode ierr;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   PetscFunctionBeginUser;
24c4762a1bSJed Brown   options->monitor          = PETSC_FALSE;
25c4762a1bSJed Brown   options->error            = PETSC_FALSE;
26c4762a1bSJed Brown   options->particlesPerCell = 1;
27c4762a1bSJed Brown   options->momentTol        = 100.0*PETSC_MACHINE_EPSILON;
28c4762a1bSJed Brown   options->omega            = 64.0;
29c4762a1bSJed Brown   options->ostep            = 100;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Harmonic Oscillator Options", "DMPLEX");CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = PetscOptionsReal("-omega", "Oscillator frequency", "ex4.c", options->omega, &options->omega, PETSC_NULL);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   PetscFunctionReturn(0);
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
42c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
43c4762a1bSJed Brown {
44c4762a1bSJed Brown   PetscErrorCode ierr;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBeginUser;
4730602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
4830602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
5130602db0SMatthew G. Knepley   ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr);
52c4762a1bSJed Brown   PetscFunctionReturn(0);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw)
56c4762a1bSJed Brown {
57c4762a1bSJed Brown   DM             dm;
58c4762a1bSJed Brown   AppCtx        *user;
59c4762a1bSJed Brown   PetscRandom    rnd;
60c4762a1bSJed Brown   PetscBool      simplex;
61c4762a1bSJed Brown   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
62c4762a1bSJed Brown   PetscInt       dim, d, cStart, cEnd, c, Np, p;
63c4762a1bSJed Brown   PetscErrorCode ierr;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   PetscFunctionBeginUser;
66c4762a1bSJed Brown   ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd);CHKERRQ(ierr);
67c4762a1bSJed Brown   ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr);
71c4762a1bSJed Brown   Np   = user->particlesPerCell;
72c4762a1bSJed Brown   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
7430602db0SMatthew G. Knepley   ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr);
77c4762a1bSJed Brown   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
78c4762a1bSJed Brown   ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
79c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
80c4762a1bSJed Brown     if (Np == 1) {
81c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
82c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
83c4762a1bSJed Brown     } else {
84c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */
85c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
86c4762a1bSJed Brown         const PetscInt n   = c*Np + p;
87c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
88c4762a1bSJed Brown 
89c4762a1bSJed Brown         for (d = 0; d < dim; ++d) {
90c4762a1bSJed Brown           ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr);
91c4762a1bSJed Brown           sum += refcoords[d];
92c4762a1bSJed Brown         }
93c4762a1bSJed Brown         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
94c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
95c4762a1bSJed Brown       }
96c4762a1bSJed Brown     }
97c4762a1bSJed Brown   }
98c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   DM             dm;
107c4762a1bSJed Brown   AppCtx        *user;
108c4762a1bSJed Brown   PetscReal     *coords;
109c4762a1bSJed Brown   PetscScalar   *initialConditions;
110c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np, p;
111c4762a1bSJed Brown   PetscErrorCode ierr;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBeginUser;
114c4762a1bSJed Brown   ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr);
115c4762a1bSJed Brown   Np   = user->particlesPerCell;
116c4762a1bSJed Brown   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr);
121c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
122c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
123c4762a1bSJed Brown       const PetscInt n = c*Np + p;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown       initialConditions[n*2+0] = DMPlex_NormD_Internal(dim, &coords[n*dim]);
126c4762a1bSJed Brown       initialConditions[n*2+1] = 0.0;
127c4762a1bSJed Brown     }
128c4762a1bSJed Brown   }
129c4762a1bSJed Brown   ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr);
130c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
13140be0ff1SMatthew G. Knepley 
13240be0ff1SMatthew G. Knepley   ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
133c4762a1bSJed Brown   PetscFunctionReturn(0);
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   PetscInt      *cellid;
139c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
140c4762a1bSJed Brown   PetscErrorCode ierr;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBeginUser;
143c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
154c4762a1bSJed Brown   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
155c4762a1bSJed Brown   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
156c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
157c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
158c4762a1bSJed Brown       const PetscInt n = c*Np + p;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown       cellid[n] = c;
161c4762a1bSJed Brown     }
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   AppCtx            *user  = (AppCtx *) ctx;
172c4762a1bSJed Brown   const PetscReal    omega = user->omega;
173c4762a1bSJed Brown   const PetscScalar *u;
174c4762a1bSJed Brown   MPI_Comm           comm;
175c4762a1bSJed Brown   PetscReal          dt;
176c4762a1bSJed Brown   PetscInt           Np, p;
177c4762a1bSJed Brown   PetscErrorCode     ierr;
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   PetscFunctionBeginUser;
180c4762a1bSJed Brown   if (step%user->ostep == 0) {
181c4762a1bSJed Brown     ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
182c4762a1bSJed Brown     if (!step) {ierr = PetscPrintf(comm, "Time     Step Part     Energy Mod Energy\n");CHKERRQ(ierr);}
183c4762a1bSJed Brown     ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr);
184c4762a1bSJed Brown     ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
185c4762a1bSJed Brown     ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
186c4762a1bSJed Brown     Np /= 2;
187c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
188c4762a1bSJed Brown       const PetscReal x  = PetscRealPart(u[p*2+0]);
189c4762a1bSJed Brown       const PetscReal v  = PetscRealPart(u[p*2+1]);
190c4762a1bSJed Brown       const PetscReal E  = 0.5*(v*v + PetscSqr(omega)*x*x);
191c4762a1bSJed Brown       const PetscReal mE = 0.5*(v*v + PetscSqr(omega)*x*x - PetscSqr(omega)*dt*x*v);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown       ierr = PetscPrintf(comm, "%.6lf %4D %4D %10.4lf %10.4lf\n", t, step, p, (double) E, (double) mE);CHKERRQ(ierr);
194c4762a1bSJed Brown     }
195c4762a1bSJed Brown     ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown   PetscFunctionReturn(0);
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   DM             dm;
203c4762a1bSJed Brown   AppCtx        *user;
204c4762a1bSJed Brown   PetscErrorCode ierr;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   PetscFunctionBeginUser;
207c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = DMGetApplicationContext(dm, (void **) &user);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = SetInitialCoordinates(dm);CHKERRQ(ierr);
210c4762a1bSJed Brown   ierr = SetInitialConditions(dm, u);CHKERRQ(ierr);
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   MPI_Comm           comm;
217c4762a1bSJed Brown   DM                 sdm;
218c4762a1bSJed Brown   AppCtx            *user;
219c4762a1bSJed Brown   const PetscScalar *u, *coords;
220c4762a1bSJed Brown   PetscScalar       *e;
221c4762a1bSJed Brown   PetscReal          t, omega;
222c4762a1bSJed Brown   PetscInt           dim, Np, p;
223c4762a1bSJed Brown   PetscErrorCode     ierr;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   PetscFunctionBeginUser;
226c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
227c4762a1bSJed Brown   ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = DMGetApplicationContext(sdm, (void **) &user);CHKERRQ(ierr);
229c4762a1bSJed Brown   omega = user->omega;
230c4762a1bSJed Brown   ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = TSGetSolveTime(ts, &t);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = VecGetArray(E, &e);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
235c4762a1bSJed Brown   Np  /= 2;
236c4762a1bSJed Brown   ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
237c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
238c4762a1bSJed Brown     const PetscReal x  = PetscRealPart(u[p*2+0]);
239c4762a1bSJed Brown     const PetscReal v  = PetscRealPart(u[p*2+1]);
240c4762a1bSJed Brown     const PetscReal x0 = DMPlex_NormD_Internal(dim, &coords[p*dim]);
241c4762a1bSJed Brown     const PetscReal ex =  x0*PetscCosReal(omega*t);
242c4762a1bSJed Brown     const PetscReal ev = -x0*omega*PetscSinReal(omega*t);
243c4762a1bSJed Brown 
244c4762a1bSJed Brown     if (user->error) {ierr = PetscPrintf(comm, "p%D error [%.2g %.2g] sol [%.6lf %.6lf] exact [%.6lf %.6lf] energy/exact energy %g / %g\n", p, (double) PetscAbsReal(x-ex), (double) PetscAbsReal(v-ev), (double) x, (double) v, (double) ex, (double) ev, 0.5*(v*v + PetscSqr(omega)*x*x), (double) 0.5*PetscSqr(omega*x0));}
245c4762a1bSJed Brown     e[p*2+0] = x - ex;
246c4762a1bSJed Brown     e[p*2+1] = v - ev;
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown   ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
24940be0ff1SMatthew G. Knepley 
250c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = VecRestoreArray(E, &e);CHKERRQ(ierr);
252c4762a1bSJed Brown   PetscFunctionReturn(0);
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
25540be0ff1SMatthew G. Knepley /*---------------------Create particle RHS Functions--------------------------*/
25640be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
25740be0ff1SMatthew G. Knepley {
25840be0ff1SMatthew G. Knepley   const PetscScalar *v;
25940be0ff1SMatthew G. Knepley   PetscScalar       *xres;
26040be0ff1SMatthew G. Knepley   PetscInt           Np, p;
26140be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
26240be0ff1SMatthew G. Knepley 
26340be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
26440be0ff1SMatthew G. Knepley   ierr = VecGetArray(Xres, &xres);CHKERRQ(ierr);
26540be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(V, &v);CHKERRQ(ierr);
26640be0ff1SMatthew G. Knepley   ierr = VecGetLocalSize(Xres, &Np);CHKERRQ(ierr);
26740be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
26840be0ff1SMatthew G. Knepley      xres[p] = v[p];
26940be0ff1SMatthew G. Knepley   }
27040be0ff1SMatthew G. Knepley   ierr = VecRestoreArrayRead(V, &v);CHKERRQ(ierr);
27140be0ff1SMatthew G. Knepley   ierr = VecRestoreArray(Xres, &xres);CHKERRQ(ierr);
27240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
27340be0ff1SMatthew G. Knepley }
27440be0ff1SMatthew G. Knepley 
27540be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
27640be0ff1SMatthew G. Knepley {
27740be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *)ctx;
27840be0ff1SMatthew G. Knepley   const PetscScalar *x;
27940be0ff1SMatthew G. Knepley   PetscScalar       *vres;
28040be0ff1SMatthew G. Knepley   PetscInt           Np, p;
28140be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
28240be0ff1SMatthew G. Knepley 
28340be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
28440be0ff1SMatthew G. Knepley   ierr = VecGetArray(Vres, &vres);CHKERRQ(ierr);
28540be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
28640be0ff1SMatthew G. Knepley   ierr = VecGetLocalSize(Vres, &Np);CHKERRQ(ierr);
28740be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
28840be0ff1SMatthew G. Knepley     vres[p] = -PetscSqr(user->omega)*x[p];
28940be0ff1SMatthew G. Knepley   }
29040be0ff1SMatthew G. Knepley   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
29140be0ff1SMatthew G. Knepley   ierr = VecRestoreArray(Vres, &vres);CHKERRQ(ierr);
29240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
29340be0ff1SMatthew G. Knepley }
29440be0ff1SMatthew G. Knepley 
29540be0ff1SMatthew G. Knepley // static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
29640be0ff1SMatthew G. Knepley // {
29740be0ff1SMatthew G. Knepley //   AppCtx            *user = (AppCtx *) ctx;
29840be0ff1SMatthew G. Knepley //   DM                 dm;
29940be0ff1SMatthew G. Knepley //   const PetscScalar *u;
30040be0ff1SMatthew G. Knepley //   PetscScalar       *r;
30140be0ff1SMatthew G. Knepley //   PetscInt           Np, p;
30240be0ff1SMatthew G. Knepley //   PetscErrorCode     ierr;
30340be0ff1SMatthew G. Knepley //
30440be0ff1SMatthew G. Knepley //   PetscFunctionBeginUser;
30540be0ff1SMatthew G. Knepley //   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
30640be0ff1SMatthew G. Knepley //   ierr = VecGetArray(R, &r);CHKERRQ(ierr);
30740be0ff1SMatthew G. Knepley //   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
30840be0ff1SMatthew G. Knepley //   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
30940be0ff1SMatthew G. Knepley //   Np  /= 2;
31040be0ff1SMatthew G. Knepley //   for (p = 0; p < Np; ++p) {
31140be0ff1SMatthew G. Knepley //     r[p*2+0] = u[p*2+1];
31240be0ff1SMatthew G. Knepley //     r[p*2+1] = -PetscSqr(user->omega)*u[p*2+0];
31340be0ff1SMatthew G. Knepley //   }
31440be0ff1SMatthew G. Knepley //   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
31540be0ff1SMatthew G. Knepley //   ierr = VecRestoreArray(R, &r);CHKERRQ(ierr);
31640be0ff1SMatthew G. Knepley //   PetscFunctionReturn(0);
31740be0ff1SMatthew G. Knepley // }
31840be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/
31940be0ff1SMatthew G. Knepley 
32040be0ff1SMatthew G. Knepley /*--------------------Define RHSFunction, RHSJacobian (Theta)-----------------*/
32140be0ff1SMatthew G. Knepley static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
32240be0ff1SMatthew G. Knepley {
32340be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
32440be0ff1SMatthew G. Knepley   DM                 dm;
32540be0ff1SMatthew G. Knepley   const PetscScalar *u;
32640be0ff1SMatthew G. Knepley   PetscScalar       *g;
32740be0ff1SMatthew G. Knepley   PetscInt           Np, p;
32840be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
32940be0ff1SMatthew G. Knepley 
33040be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
33140be0ff1SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
33240be0ff1SMatthew G. Knepley   ierr = VecGetArray(G, &g);CHKERRQ(ierr);
33340be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
33440be0ff1SMatthew G. Knepley   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
33540be0ff1SMatthew G. Knepley   Np  /= 2;
33640be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
33740be0ff1SMatthew G. Knepley     g[p*2+0] = u[p*2+1];
33840be0ff1SMatthew G. Knepley     g[p*2+1] = -PetscSqr(user->omega)*u[p*2+0];
33940be0ff1SMatthew G. Knepley   }
34040be0ff1SMatthew G. Knepley   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
34140be0ff1SMatthew G. Knepley   ierr = VecRestoreArray(G, &g);CHKERRQ(ierr);
34240be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
34340be0ff1SMatthew G. Knepley }
34440be0ff1SMatthew G. Knepley 
34540be0ff1SMatthew G. Knepley /*Ji = dFi/dxj
34640be0ff1SMatthew G. Knepley J= (0    1)
34740be0ff1SMatthew G. Knepley    (-w^2 0)
34840be0ff1SMatthew G. Knepley */
34940be0ff1SMatthew G. Knepley static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U , Mat J, Mat P, void *ctx)
35040be0ff1SMatthew G. Knepley {
35140be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
35240be0ff1SMatthew G. Knepley   PetscInt           Np = user->dim * user->particlesPerCell;
35340be0ff1SMatthew G. Knepley   PetscInt           i, m, n;
35440be0ff1SMatthew G. Knepley   const PetscScalar *u;
35540be0ff1SMatthew G. Knepley   PetscScalar        vals[4] = {0., 1., -PetscSqr(user->omega), 0.};
35640be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
35740be0ff1SMatthew G. Knepley 
358*362febeeSStefano Zampini   PetscFunctionBeginUser;
35940be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
36040be0ff1SMatthew G. Knepley   //ierr = PetscPrintf(PETSC_COMM_WORLD, "# Particles (Np) = %d \n" , Np);CHKERRQ(ierr);
36140be0ff1SMatthew G. Knepley 
36240be0ff1SMatthew G. Knepley   ierr = MatGetOwnershipRange(J, &m, &n);CHKERRQ(ierr);
36340be0ff1SMatthew G. Knepley   for (i = 0; i < Np; ++i) {
36440be0ff1SMatthew G. Knepley     const PetscInt rows[2] = {2*i, 2*i+1};
36540be0ff1SMatthew G. Knepley     ierr = MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES);CHKERRQ(ierr);
36640be0ff1SMatthew G. Knepley   }
36740be0ff1SMatthew G. Knepley   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36840be0ff1SMatthew G. Knepley   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36940be0ff1SMatthew G. Knepley   //ierr = MatView(J,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
37040be0ff1SMatthew G. Knepley 
37140be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
37240be0ff1SMatthew G. Knepley 
37340be0ff1SMatthew G. Knepley }
37440be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/
37540be0ff1SMatthew G. Knepley 
37640be0ff1SMatthew G. Knepley /*----------------Define S, F, G Functions (Discrete Gradients)---------------*/
37740be0ff1SMatthew G. Knepley /*
37840be0ff1SMatthew G. Knepley   u_t = S * gradF
37940be0ff1SMatthew G. Knepley     --or--
38040be0ff1SMatthew G. Knepley   u_t = S * G
38140be0ff1SMatthew G. Knepley 
38240be0ff1SMatthew G. Knepley   + Sfunc - constructor for the S matrix from the formulation
38340be0ff1SMatthew G. Knepley   . Ffunc - functional F from the formulation
38440be0ff1SMatthew G. Knepley   - Gfunc - constructor for the gradient of F from the formulation
38540be0ff1SMatthew G. Knepley */
38640be0ff1SMatthew G. Knepley 
38740be0ff1SMatthew G. Knepley PetscErrorCode Sfunc(TS ts, PetscReal t, Vec U, Mat S, void *ctx)
38840be0ff1SMatthew G. Knepley {
38940be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
39040be0ff1SMatthew G. Knepley   PetscInt           Np = user->dim * user->particlesPerCell;
39140be0ff1SMatthew G. Knepley   PetscInt           i, m, n;
39240be0ff1SMatthew G. Knepley   const PetscScalar *u;
39340be0ff1SMatthew G. Knepley   PetscScalar        vals[4] = {0., 1., -1, 0.};
39440be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
39540be0ff1SMatthew G. Knepley 
39640be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
39740be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
39840be0ff1SMatthew G. Knepley   ierr = MatGetOwnershipRange(S, &m, &n);CHKERRQ(ierr);
39940be0ff1SMatthew G. Knepley   for (i = 0; i < Np; ++i) {
40040be0ff1SMatthew G. Knepley     const PetscInt rows[2] = {2*i, 2*i+1};
40140be0ff1SMatthew G. Knepley     ierr = MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES);CHKERRQ(ierr);
40240be0ff1SMatthew G. Knepley   }
40340be0ff1SMatthew G. Knepley   ierr = MatAssemblyBegin(S,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40440be0ff1SMatthew G. Knepley   ierr = MatAssemblyEnd(S,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
40540be0ff1SMatthew G. Knepley 
40640be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
40740be0ff1SMatthew G. Knepley }
40840be0ff1SMatthew G. Knepley 
40940be0ff1SMatthew G. Knepley PetscErrorCode Ffunc(TS ts, PetscReal t, Vec U, PetscScalar *F, void *ctx)
41040be0ff1SMatthew G. Knepley {
41140be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
41240be0ff1SMatthew G. Knepley   DM                 dm;
41340be0ff1SMatthew G. Knepley   const PetscScalar *u;
41440be0ff1SMatthew G. Knepley   PetscInt           Np = user->dim * user->particlesPerCell;
41540be0ff1SMatthew G. Knepley   PetscInt           p;
41640be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
41740be0ff1SMatthew G. Knepley 
41840be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
41940be0ff1SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
42040be0ff1SMatthew G. Knepley   //Define F
42140be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
42240be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
42340be0ff1SMatthew G. Knepley     *F += (1/2)*PetscSqr(user->omega)*PetscSqr(u[p*2+0]) + (1/2)*PetscSqr(u[p*2+1]);
42440be0ff1SMatthew G. Knepley   }
42540be0ff1SMatthew G. Knepley   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
42640be0ff1SMatthew G. Knepley 
42740be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
42840be0ff1SMatthew G. Knepley }
42940be0ff1SMatthew G. Knepley 
43040be0ff1SMatthew G. Knepley PetscErrorCode gradFfunc(TS ts, PetscReal t, Vec U, Vec gradF, void *ctx)
43140be0ff1SMatthew G. Knepley {
43240be0ff1SMatthew G. Knepley   AppCtx            *user = (AppCtx *) ctx;
43340be0ff1SMatthew G. Knepley   DM                 dm;
43440be0ff1SMatthew G. Knepley   const PetscScalar *u;
43540be0ff1SMatthew G. Knepley   PetscScalar       *g;
43640be0ff1SMatthew G. Knepley   PetscInt           Np = user->dim * user->particlesPerCell;
43740be0ff1SMatthew G. Knepley   PetscInt           p;
43840be0ff1SMatthew G. Knepley   PetscErrorCode     ierr;
43940be0ff1SMatthew G. Knepley 
44040be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
44140be0ff1SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
44240be0ff1SMatthew G. Knepley   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
44340be0ff1SMatthew G. Knepley 
44440be0ff1SMatthew G. Knepley   //Define gradF
44540be0ff1SMatthew G. Knepley   ierr = VecGetArray(gradF, &g);CHKERRQ(ierr);
44640be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
44740be0ff1SMatthew G. Knepley     g[p*2+0] = PetscSqr(user->omega)*u[p*2+0]; /*dF/dx*/
44840be0ff1SMatthew G. Knepley     g[p*2+1] = u[p*2+1]; /*dF/dv*/
44940be0ff1SMatthew G. Knepley   }
45040be0ff1SMatthew G. Knepley   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
45140be0ff1SMatthew G. Knepley   ierr = VecRestoreArray(gradF, &g);CHKERRQ(ierr);
45240be0ff1SMatthew G. Knepley 
45340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
45440be0ff1SMatthew G. Knepley }
45540be0ff1SMatthew G. Knepley /*----------------------------------------------------------------------------*/
45640be0ff1SMatthew G. Knepley 
457c4762a1bSJed Brown int main(int argc,char **argv)
458c4762a1bSJed Brown {
459c4762a1bSJed Brown   TS             ts;     /* nonlinear solver                 */
460c4762a1bSJed Brown   DM             dm, sw; /* Mesh and particle managers       */
461c4762a1bSJed Brown   Vec            u;      /* swarm vector                     */
46240be0ff1SMatthew G. Knepley   PetscInt       n;      /* swarm vector size                */
463c4762a1bSJed Brown   IS             is1, is2;
464c4762a1bSJed Brown   MPI_Comm       comm;
46540be0ff1SMatthew G. Knepley   Mat            J;      /* Jacobian matrix                  */
466c4762a1bSJed Brown   AppCtx         user;
467c4762a1bSJed Brown   PetscErrorCode ierr;
468c4762a1bSJed Brown 
46940be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47040be0ff1SMatthew G. Knepley      Initialize program
47140be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
472c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
473c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
474c4762a1bSJed Brown   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
475c4762a1bSJed Brown 
47640be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47740be0ff1SMatthew G. Knepley      Create Particle-Mesh
47840be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
479c4762a1bSJed Brown   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
480c4762a1bSJed Brown   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
481c4762a1bSJed Brown   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
482c4762a1bSJed Brown 
48340be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48440be0ff1SMatthew G. Knepley      Setup timestepping solver
48540be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
486c4762a1bSJed Brown   ierr = TSCreate(comm, &ts);CHKERRQ(ierr);
48740be0ff1SMatthew G. Knepley   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
48840be0ff1SMatthew G. Knepley 
489c4762a1bSJed Brown   ierr = TSSetDM(ts, sw);CHKERRQ(ierr);
490c4762a1bSJed Brown   ierr = TSSetMaxTime(ts, 0.1);CHKERRQ(ierr);
491c4762a1bSJed Brown   ierr = TSSetTimeStep(ts, 0.00001);CHKERRQ(ierr);
492c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts, 100);CHKERRQ(ierr);
493c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
494c4762a1bSJed Brown   if (user.monitor) {ierr = TSMonitorSet(ts, Monitor, &user, NULL);CHKERRQ(ierr);}
495c4762a1bSJed Brown 
49640be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49740be0ff1SMatthew G. Knepley      Prepare to solve
49840be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
499c4762a1bSJed Brown   ierr = DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr);
500c4762a1bSJed Brown   ierr = VecGetLocalSize(u, &n);CHKERRQ(ierr);
50140be0ff1SMatthew G. Knepley   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
50240be0ff1SMatthew G. Knepley   ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr);
50340be0ff1SMatthew G. Knepley   ierr = TSSetComputeExactError(ts, ComputeError);CHKERRQ(ierr);
50440be0ff1SMatthew G. Knepley 
50540be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50640be0ff1SMatthew G. Knepley      Define function F(U, Udot , x , t) = G(U, x, t)
50740be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50840be0ff1SMatthew G. Knepley 
50940be0ff1SMatthew G. Knepley   /* - - - - - - - Basic Symplectic - - - - - - - - - - - - - - - - - - - - - -*/
510c4762a1bSJed Brown   ierr = ISCreateStride(comm, n/2, 0, 2, &is1);CHKERRQ(ierr);
511c4762a1bSJed Brown   ierr = ISCreateStride(comm, n/2, 1, 2, &is2);CHKERRQ(ierr);
512c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts, "position", is1);CHKERRQ(ierr);
513c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts, "momentum", is2);CHKERRQ(ierr);
514c4762a1bSJed Brown   ierr = ISDestroy(&is1);CHKERRQ(ierr);
515c4762a1bSJed Brown   ierr = ISDestroy(&is2);CHKERRQ(ierr);
516c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunction1, &user);CHKERRQ(ierr);
517c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunction2, &user);CHKERRQ(ierr);
518c4762a1bSJed Brown 
51940be0ff1SMatthew G. Knepley   /* - - - - - - - Theta (Implicit Midpoint) - - - - - - - - - - - - - - - - -*/
52040be0ff1SMatthew G. Knepley   ierr = TSSetRHSFunction(ts, NULL, RHSFunction, &user);CHKERRQ(ierr);
52140be0ff1SMatthew G. Knepley 
52240be0ff1SMatthew G. Knepley   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
52340be0ff1SMatthew G. Knepley   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
52440be0ff1SMatthew G. Knepley   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
52540be0ff1SMatthew G. Knepley   ierr = MatSetUp(J);CHKERRQ(ierr);
52640be0ff1SMatthew G. Knepley   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
52740be0ff1SMatthew G. Knepley   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
52840be0ff1SMatthew G. Knepley   ierr = PetscPrintf(comm, "n = %d\n",n);CHKERRQ(ierr);//Check number of particles
52940be0ff1SMatthew G. Knepley   ierr = TSSetRHSJacobian(ts,J,J,RHSJacobian,&user);CHKERRQ(ierr);
53040be0ff1SMatthew G. Knepley 
53140be0ff1SMatthew G. Knepley   /* - - - - - - - Discrete Gradients - - - - - - - - - - - - - - - - - - - - */
53240be0ff1SMatthew G. Knepley   ierr = TSDiscGradSetFormulation(ts, Sfunc, Ffunc, gradFfunc, &user);CHKERRQ(ierr);
53340be0ff1SMatthew G. Knepley 
53440be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53540be0ff1SMatthew G. Knepley      Solve
53640be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53740be0ff1SMatthew G. Knepley 
538c4762a1bSJed Brown   ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr);
539c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
540c4762a1bSJed Brown 
54140be0ff1SMatthew G. Knepley   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54240be0ff1SMatthew G. Knepley      Clean up workspace
54340be0ff1SMatthew G. Knepley      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54440be0ff1SMatthew G. Knepley   ierr = MatDestroy(&J);CHKERRQ(ierr);
54540be0ff1SMatthew G. Knepley   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr);
546c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
547c4762a1bSJed Brown   ierr = DMDestroy(&sw);CHKERRQ(ierr);
548c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
549c4762a1bSJed Brown   ierr = PetscFinalize();
550c4762a1bSJed Brown   return ierr;
551c4762a1bSJed Brown }
552c4762a1bSJed Brown 
553c4762a1bSJed Brown /*TEST
554c4762a1bSJed Brown 
555c4762a1bSJed Brown    build:
556c4762a1bSJed Brown      requires: triangle !single !complex
557c4762a1bSJed Brown    test:
558c4762a1bSJed Brown      suffix: 1
55940be0ff1SMatthew G. Knepley      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 1 -ts_convergence_estimate -convest_num_refine 2 -dm_view -sw_view -monitor -output_step 50 -error
56040be0ff1SMatthew G. Knepley 
561c4762a1bSJed Brown    test:
562c4762a1bSJed Brown      suffix: 2
56340be0ff1SMatthew G. Knepley      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 2 -ts_convergence_estimate -convest_num_refine 2 -dm_view -sw_view -monitor -output_step 50 -error
56440be0ff1SMatthew G. Knepley 
565c4762a1bSJed Brown    test:
566c4762a1bSJed Brown      suffix: 3
56740be0ff1SMatthew G. Knepley      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -sw_view -monitor -output_step 50 -error
56840be0ff1SMatthew G. Knepley 
569c4762a1bSJed Brown    test:
570c4762a1bSJed Brown      suffix: 4
57140be0ff1SMatthew G. Knepley      args: -dm_plex_box_faces 1,1 -ts_type basicsymplectic -ts_basicsymplectic_type 4 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -sw_view -monitor -output_step 50 -error
57240be0ff1SMatthew G. Knepley 
57340be0ff1SMatthew G. Knepley    test:
57440be0ff1SMatthew G. Knepley      suffix: 5
57540be0ff1SMatthew G. Knepley      args: -dm_plex_box_faces 1,1 -ts_type theta -ts_theta_theta 0.5 -ts_convergence_estimate -convest_num_refine 2 -dm_view -sw_view -monitor -output_step 50 -error
57640be0ff1SMatthew G. Knepley 
57740be0ff1SMatthew G. Knepley    test:
57840be0ff1SMatthew G. Knepley      suffix: 6
57940be0ff1SMatthew G. Knepley      args: -dm_plex_box_faces 1,1 -ts_type discgrad -monitor -output_step 50 -error -ts_convergence_estimate -convest_num_refine 2
580c4762a1bSJed Brown 
581c4762a1bSJed Brown TEST*/
582