xref: /petsc/src/dm/impls/swarm/tests/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Example of simple hamiltonian system (harmonic oscillator) with particles and a basic symplectic integrator\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmplex.h>
4*c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>  /* For norm */
5*c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> /* For CoordinatesRefToReal() */
6*c4762a1bSJed Brown #include <petscdmswarm.h>
7*c4762a1bSJed Brown #include <petscts.h>
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown typedef struct {
10*c4762a1bSJed Brown   PetscInt  dim;                          /* The topological mesh dimension */
11*c4762a1bSJed Brown   PetscBool simplex;                      /* Flag for simplices or tensor cells */
12*c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */
13*c4762a1bSJed Brown   PetscReal omega;                        /* Oscillation frequency omega */
14*c4762a1bSJed Brown   PetscInt  particlesPerCell;             /* The number of partices per cell */
15*c4762a1bSJed Brown   PetscReal momentTol;                    /* Tolerance for checking moment conservation */
16*c4762a1bSJed Brown   PetscBool monitor;                      /* Flag for using the TS monitor */
17*c4762a1bSJed Brown   PetscBool error;                        /* Flag for printing the error */
18*c4762a1bSJed Brown   PetscInt  ostep;                        /* print the energy at each ostep time steps */
19*c4762a1bSJed Brown } AppCtx;
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22*c4762a1bSJed Brown {
23*c4762a1bSJed Brown   PetscErrorCode ierr;
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   PetscFunctionBeginUser;
26*c4762a1bSJed Brown   options->dim              = 2;
27*c4762a1bSJed Brown   options->simplex          = PETSC_TRUE;
28*c4762a1bSJed Brown   options->monitor          = PETSC_FALSE;
29*c4762a1bSJed Brown   options->error            = PETSC_FALSE;
30*c4762a1bSJed Brown   options->particlesPerCell = 1;
31*c4762a1bSJed Brown   options->momentTol        = 100.0*PETSC_MACHINE_EPSILON;
32*c4762a1bSJed Brown   options->omega            = 64.0;
33*c4762a1bSJed Brown   options->ostep            = 100;
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   ierr = PetscStrcpy(options->filename, "");CHKERRQ(ierr);
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Harmonic Oscillator Options", "DMPLEX");CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = PetscOptionsBool("-monitor", "Flag to use the TS monitor", "ex4.c", options->monitor, &options->monitor, NULL);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "The flag for simplices or tensor cells", "ex4.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscOptionsString("-mesh", "Name of the mesh filename if any", "ex4.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = PetscOptionsInt("-particles_per_cell", "Number of particles per cell", "ex4.c", options->particlesPerCell, &options->particlesPerCell, NULL);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = PetscOptionsReal("-omega", "Oscillator frequency", "ex4.c", options->omega, &options->omega, PETSC_NULL);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   PetscFunctionReturn(0);
49*c4762a1bSJed Brown }
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
52*c4762a1bSJed Brown {
53*c4762a1bSJed Brown   PetscBool      flg;
54*c4762a1bSJed Brown   PetscErrorCode ierr;
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   PetscFunctionBeginUser;
57*c4762a1bSJed Brown   ierr = PetscStrcmp(user->filename, "", &flg);CHKERRQ(ierr);
58*c4762a1bSJed Brown   if (flg) {
59*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
60*c4762a1bSJed Brown   } else {
61*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr);
62*c4762a1bSJed Brown     ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr);
63*c4762a1bSJed Brown   }
64*c4762a1bSJed Brown   {
65*c4762a1bSJed Brown     DM distributedMesh = NULL;
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
68*c4762a1bSJed Brown     if (distributedMesh) {
69*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
70*c4762a1bSJed Brown       *dm  = distributedMesh;
71*c4762a1bSJed Brown     }
72*c4762a1bSJed Brown   }
73*c4762a1bSJed Brown   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */
74*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
77*c4762a1bSJed Brown   PetscFunctionReturn(0);
78*c4762a1bSJed Brown }
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown static PetscErrorCode SetInitialCoordinates(DM dmSw)
81*c4762a1bSJed Brown {
82*c4762a1bSJed Brown   DM             dm;
83*c4762a1bSJed Brown   AppCtx        *user;
84*c4762a1bSJed Brown   PetscRandom    rnd;
85*c4762a1bSJed Brown   PetscBool      simplex;
86*c4762a1bSJed Brown   PetscReal     *centroid, *coords, *xi0, *v0, *J, *invJ, detJ;
87*c4762a1bSJed Brown   PetscInt       dim, d, cStart, cEnd, c, Np, p;
88*c4762a1bSJed Brown   PetscErrorCode ierr;
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   PetscFunctionBeginUser;
91*c4762a1bSJed Brown   ierr = PetscRandomCreate(PetscObjectComm((PetscObject) dmSw), &rnd);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = PetscRandomSetInterval(rnd, -1.0, 1.0);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr);
96*c4762a1bSJed Brown   simplex = user->simplex;
97*c4762a1bSJed Brown   Np   = user->particlesPerCell;
98*c4762a1bSJed Brown   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = PetscMalloc5(dim, &centroid, dim, &xi0, dim, &v0, dim*dim, &J, dim*dim, &invJ);CHKERRQ(ierr);
102*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) xi0[d] = -1.0;
103*c4762a1bSJed Brown   ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
104*c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
105*c4762a1bSJed Brown     if (Np == 1) {
106*c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL);CHKERRQ(ierr);
107*c4762a1bSJed Brown       for (d = 0; d < dim; ++d) coords[c*dim+d] = centroid[d];
108*c4762a1bSJed Brown     } else {
109*c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); /* affine */
110*c4762a1bSJed Brown       for (p = 0; p < Np; ++p) {
111*c4762a1bSJed Brown         const PetscInt n   = c*Np + p;
112*c4762a1bSJed Brown         PetscReal      sum = 0.0, refcoords[3];
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown         for (d = 0; d < dim; ++d) {
115*c4762a1bSJed Brown           ierr = PetscRandomGetValueReal(rnd, &refcoords[d]);CHKERRQ(ierr);
116*c4762a1bSJed Brown           sum += refcoords[d];
117*c4762a1bSJed Brown         }
118*c4762a1bSJed Brown         if (simplex && sum > 0.0) for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim)*sum;
119*c4762a1bSJed Brown         CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n*dim]);
120*c4762a1bSJed Brown       }
121*c4762a1bSJed Brown     }
122*c4762a1bSJed Brown   }
123*c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = PetscFree5(centroid, xi0, v0, J, invJ);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
126*c4762a1bSJed Brown   PetscFunctionReturn(0);
127*c4762a1bSJed Brown }
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown static PetscErrorCode SetInitialConditions(DM dmSw, Vec u)
130*c4762a1bSJed Brown {
131*c4762a1bSJed Brown   DM             dm;
132*c4762a1bSJed Brown   AppCtx        *user;
133*c4762a1bSJed Brown   PetscReal     *coords;
134*c4762a1bSJed Brown   PetscScalar   *initialConditions;
135*c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np, p;
136*c4762a1bSJed Brown   PetscErrorCode ierr;
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   PetscFunctionBeginUser;
139*c4762a1bSJed Brown   ierr = DMGetApplicationContext(dmSw, (void **) &user);CHKERRQ(ierr);
140*c4762a1bSJed Brown   Np   = user->particlesPerCell;
141*c4762a1bSJed Brown   ierr = DMSwarmGetCellDM(dmSw, &dm);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = DMSwarmGetField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = VecGetArray(u, &initialConditions);CHKERRQ(ierr);
146*c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
147*c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
148*c4762a1bSJed Brown       const PetscInt n = c*Np + p;
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown       initialConditions[n*2+0] = DMPlex_NormD_Internal(dim, &coords[n*dim]);
151*c4762a1bSJed Brown       initialConditions[n*2+1] = 0.0;
152*c4762a1bSJed Brown     }
153*c4762a1bSJed Brown   }
154*c4762a1bSJed Brown   ierr = VecRestoreArray(u, &initialConditions);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = DMSwarmRestoreField(dmSw, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
156*c4762a1bSJed Brown   PetscFunctionReturn(0);
157*c4762a1bSJed Brown }
158*c4762a1bSJed Brown 
159*c4762a1bSJed Brown static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user)
160*c4762a1bSJed Brown {
161*c4762a1bSJed Brown   PetscInt      *cellid;
162*c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c, Np = user->particlesPerCell, p;
163*c4762a1bSJed Brown   PetscErrorCode ierr;
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown   PetscFunctionBeginUser;
166*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = DMCreate(PetscObjectComm((PetscObject) dm), sw);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr);
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown   ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(*sw, "kinematics", 2, PETSC_REAL);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Np, 0);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = DMSetFromOptions(*sw);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
179*c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
180*c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
181*c4762a1bSJed Brown       const PetscInt n = c*Np + p;
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown       cellid[n] = c;
184*c4762a1bSJed Brown     }
185*c4762a1bSJed Brown   }
186*c4762a1bSJed Brown   ierr = DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **) &cellid);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *sw, "Particles");CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = DMViewFromOptions(*sw, NULL, "-sw_view");CHKERRQ(ierr);
189*c4762a1bSJed Brown   PetscFunctionReturn(0);
190*c4762a1bSJed Brown }
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown /* Create particle RHS Functions */
193*c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
194*c4762a1bSJed Brown {
195*c4762a1bSJed Brown   const PetscScalar *v;
196*c4762a1bSJed Brown   PetscScalar       *xres;
197*c4762a1bSJed Brown   PetscInt           Np, p;
198*c4762a1bSJed Brown   PetscErrorCode     ierr;
199*c4762a1bSJed Brown 
200*c4762a1bSJed Brown   PetscFunctionBeginUser;
201*c4762a1bSJed Brown   ierr = VecGetArray(Xres, &xres);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = VecGetArrayRead(V, &v);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = VecGetLocalSize(Xres, &Np);CHKERRQ(ierr);
204*c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
205*c4762a1bSJed Brown      xres[p] = v[p];
206*c4762a1bSJed Brown   }
207*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(V, &v);CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = VecRestoreArray(Xres, &xres);CHKERRQ(ierr);
209*c4762a1bSJed Brown   PetscFunctionReturn(0);
210*c4762a1bSJed Brown }
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
213*c4762a1bSJed Brown {
214*c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ctx;
215*c4762a1bSJed Brown   const PetscScalar *x;
216*c4762a1bSJed Brown   PetscScalar       *vres;
217*c4762a1bSJed Brown   PetscInt           Np, p;
218*c4762a1bSJed Brown   PetscErrorCode     ierr;
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown   PetscFunctionBeginUser;
221*c4762a1bSJed Brown   ierr = VecGetArray(Vres, &vres);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = VecGetLocalSize(Vres, &Np);CHKERRQ(ierr);
224*c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
225*c4762a1bSJed Brown     vres[p] = -PetscSqr(user->omega)*x[p];
226*c4762a1bSJed Brown   }
227*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = VecRestoreArray(Vres, &vres);CHKERRQ(ierr);
229*c4762a1bSJed Brown   PetscFunctionReturn(0);
230*c4762a1bSJed Brown }
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown static PetscErrorCode RHSFunctionParticles(TS ts, PetscReal t, Vec U, Vec R, void *ctx)
233*c4762a1bSJed Brown {
234*c4762a1bSJed Brown   AppCtx            *user = (AppCtx *) ctx;
235*c4762a1bSJed Brown   DM                 dm;
236*c4762a1bSJed Brown   const PetscScalar *u;
237*c4762a1bSJed Brown   PetscScalar       *r;
238*c4762a1bSJed Brown   PetscInt           Np, p;
239*c4762a1bSJed Brown   PetscErrorCode     ierr;
240*c4762a1bSJed Brown 
241*c4762a1bSJed Brown   PetscFunctionBeginUser;
242*c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = VecGetArray(R, &r);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
246*c4762a1bSJed Brown   Np  /= 2;
247*c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
248*c4762a1bSJed Brown     r[p*2+0] = u[p*2+1];
249*c4762a1bSJed Brown     r[p*2+1] = -PetscSqr(user->omega)*u[p*2+0];
250*c4762a1bSJed Brown   }
251*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
252*c4762a1bSJed Brown   ierr = VecRestoreArray(R, &r);CHKERRQ(ierr);
253*c4762a1bSJed Brown   PetscFunctionReturn(0);
254*c4762a1bSJed Brown }
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
257*c4762a1bSJed Brown {
258*c4762a1bSJed Brown   AppCtx            *user  = (AppCtx *) ctx;
259*c4762a1bSJed Brown   const PetscReal    omega = user->omega;
260*c4762a1bSJed Brown   const PetscScalar *u;
261*c4762a1bSJed Brown   MPI_Comm           comm;
262*c4762a1bSJed Brown   PetscReal          dt;
263*c4762a1bSJed Brown   PetscInt           Np, p;
264*c4762a1bSJed Brown   PetscErrorCode     ierr;
265*c4762a1bSJed Brown 
266*c4762a1bSJed Brown   PetscFunctionBeginUser;
267*c4762a1bSJed Brown   if (step%user->ostep == 0) {
268*c4762a1bSJed Brown     ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
269*c4762a1bSJed Brown     if (!step) {ierr = PetscPrintf(comm, "Time     Step Part     Energy Mod Energy\n");CHKERRQ(ierr);}
270*c4762a1bSJed Brown     ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr);
271*c4762a1bSJed Brown     ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
272*c4762a1bSJed Brown     ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
273*c4762a1bSJed Brown     Np /= 2;
274*c4762a1bSJed Brown     for (p = 0; p < Np; ++p) {
275*c4762a1bSJed Brown       const PetscReal x  = PetscRealPart(u[p*2+0]);
276*c4762a1bSJed Brown       const PetscReal v  = PetscRealPart(u[p*2+1]);
277*c4762a1bSJed Brown       const PetscReal E  = 0.5*(v*v + PetscSqr(omega)*x*x);
278*c4762a1bSJed Brown       const PetscReal mE = 0.5*(v*v + PetscSqr(omega)*x*x - PetscSqr(omega)*dt*x*v);
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown       ierr = PetscPrintf(comm, "%.6lf %4D %4D %10.4lf %10.4lf\n", t, step, p, (double) E, (double) mE);CHKERRQ(ierr);
281*c4762a1bSJed Brown     }
282*c4762a1bSJed Brown     ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
283*c4762a1bSJed Brown   }
284*c4762a1bSJed Brown   PetscFunctionReturn(0);
285*c4762a1bSJed Brown }
286*c4762a1bSJed Brown 
287*c4762a1bSJed Brown static PetscErrorCode InitializeSolve(TS ts, Vec u)
288*c4762a1bSJed Brown {
289*c4762a1bSJed Brown   DM             dm;
290*c4762a1bSJed Brown   AppCtx        *user;
291*c4762a1bSJed Brown   PetscErrorCode ierr;
292*c4762a1bSJed Brown 
293*c4762a1bSJed Brown   PetscFunctionBeginUser;
294*c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
295*c4762a1bSJed Brown   ierr = DMGetApplicationContext(dm, (void **) &user);CHKERRQ(ierr);
296*c4762a1bSJed Brown   ierr = SetInitialCoordinates(dm);CHKERRQ(ierr);
297*c4762a1bSJed Brown   ierr = SetInitialConditions(dm, u);CHKERRQ(ierr);
298*c4762a1bSJed Brown   PetscFunctionReturn(0);
299*c4762a1bSJed Brown }
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
302*c4762a1bSJed Brown {
303*c4762a1bSJed Brown   MPI_Comm           comm;
304*c4762a1bSJed Brown   DM                 sdm;
305*c4762a1bSJed Brown   AppCtx            *user;
306*c4762a1bSJed Brown   const PetscScalar *u, *coords;
307*c4762a1bSJed Brown   PetscScalar       *e;
308*c4762a1bSJed Brown   PetscReal          t, omega;
309*c4762a1bSJed Brown   PetscInt           dim, Np, p;
310*c4762a1bSJed Brown   PetscErrorCode     ierr;
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown   PetscFunctionBeginUser;
313*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) ts, &comm);CHKERRQ(ierr);
314*c4762a1bSJed Brown   ierr = TSGetDM(ts, &sdm);CHKERRQ(ierr);
315*c4762a1bSJed Brown   ierr = DMGetApplicationContext(sdm, (void **) &user);CHKERRQ(ierr);
316*c4762a1bSJed Brown   omega = user->omega;
317*c4762a1bSJed Brown   ierr = DMGetDimension(sdm, &dim);CHKERRQ(ierr);
318*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts, &t);CHKERRQ(ierr);
319*c4762a1bSJed Brown   ierr = VecGetArray(E, &e);CHKERRQ(ierr);
320*c4762a1bSJed Brown   ierr = VecGetArrayRead(U, &u);CHKERRQ(ierr);
321*c4762a1bSJed Brown   ierr = VecGetLocalSize(U, &Np);CHKERRQ(ierr);
322*c4762a1bSJed Brown   Np  /= 2;
323*c4762a1bSJed Brown   ierr = DMSwarmGetField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
324*c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
325*c4762a1bSJed Brown     const PetscReal x  = PetscRealPart(u[p*2+0]);
326*c4762a1bSJed Brown     const PetscReal v  = PetscRealPart(u[p*2+1]);
327*c4762a1bSJed Brown     const PetscReal x0 = DMPlex_NormD_Internal(dim, &coords[p*dim]);
328*c4762a1bSJed Brown     const PetscReal ex =  x0*PetscCosReal(omega*t);
329*c4762a1bSJed Brown     const PetscReal ev = -x0*omega*PetscSinReal(omega*t);
330*c4762a1bSJed Brown 
331*c4762a1bSJed 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));}
332*c4762a1bSJed Brown     e[p*2+0] = x - ex;
333*c4762a1bSJed Brown     e[p*2+1] = v - ev;
334*c4762a1bSJed Brown   }
335*c4762a1bSJed Brown   ierr = DMSwarmRestoreField(sdm, DMSwarmPICField_coor, NULL, NULL, (void **) &coords);CHKERRQ(ierr);
336*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U, &u);CHKERRQ(ierr);
337*c4762a1bSJed Brown   ierr = VecRestoreArray(E, &e);CHKERRQ(ierr);
338*c4762a1bSJed Brown   PetscFunctionReturn(0);
339*c4762a1bSJed Brown }
340*c4762a1bSJed Brown 
341*c4762a1bSJed Brown int main(int argc,char **argv)
342*c4762a1bSJed Brown {
343*c4762a1bSJed Brown   TS             ts;     /* nonlinear solver */
344*c4762a1bSJed Brown   DM             dm, sw; /* Mesh and particle managers */
345*c4762a1bSJed Brown   Vec            u;      /* swarm vector */
346*c4762a1bSJed Brown   IS             is1, is2;
347*c4762a1bSJed Brown   PetscInt       n;
348*c4762a1bSJed Brown   MPI_Comm       comm;
349*c4762a1bSJed Brown   AppCtx         user;
350*c4762a1bSJed Brown   PetscErrorCode ierr;
351*c4762a1bSJed Brown 
352*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
353*c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
354*c4762a1bSJed Brown   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
355*c4762a1bSJed Brown 
356*c4762a1bSJed Brown   ierr = CreateMesh(comm, &dm, &user);CHKERRQ(ierr);
357*c4762a1bSJed Brown   ierr = CreateParticles(dm, &sw, &user);CHKERRQ(ierr);
358*c4762a1bSJed Brown   ierr = DMSetApplicationContext(sw, &user);CHKERRQ(ierr);
359*c4762a1bSJed Brown 
360*c4762a1bSJed Brown   ierr = TSCreate(comm, &ts);CHKERRQ(ierr);
361*c4762a1bSJed Brown   ierr = TSSetType(ts, TSBASICSYMPLECTIC);CHKERRQ(ierr);
362*c4762a1bSJed Brown   ierr = TSSetDM(ts, sw);CHKERRQ(ierr);
363*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts, 0.1);CHKERRQ(ierr);
364*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts, 0.00001);CHKERRQ(ierr);
365*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts, 100);CHKERRQ(ierr);
366*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
367*c4762a1bSJed Brown   if (user.monitor) {ierr = TSMonitorSet(ts, Monitor, &user, NULL);CHKERRQ(ierr);}
368*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts, NULL, RHSFunctionParticles, &user);CHKERRQ(ierr);
369*c4762a1bSJed Brown 
370*c4762a1bSJed Brown   ierr = DMSwarmCreateGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr);
371*c4762a1bSJed Brown   ierr = VecGetLocalSize(u, &n);CHKERRQ(ierr);
372*c4762a1bSJed Brown   ierr = ISCreateStride(comm, n/2, 0, 2, &is1);CHKERRQ(ierr);
373*c4762a1bSJed Brown   ierr = ISCreateStride(comm, n/2, 1, 2, &is2);CHKERRQ(ierr);
374*c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts, "position", is1);CHKERRQ(ierr);
375*c4762a1bSJed Brown   ierr = TSRHSSplitSetIS(ts, "momentum", is2);CHKERRQ(ierr);
376*c4762a1bSJed Brown   ierr = ISDestroy(&is1);CHKERRQ(ierr);
377*c4762a1bSJed Brown   ierr = ISDestroy(&is2);CHKERRQ(ierr);
378*c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunction1, &user);CHKERRQ(ierr);
379*c4762a1bSJed Brown   ierr = TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunction2, &user);CHKERRQ(ierr);
380*c4762a1bSJed Brown 
381*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
382*c4762a1bSJed Brown   ierr = TSSetComputeInitialCondition(ts, InitializeSolve);CHKERRQ(ierr);
383*c4762a1bSJed Brown   ierr = TSSetComputeExactError(ts, ComputeError);CHKERRQ(ierr);
384*c4762a1bSJed Brown   ierr = TSComputeInitialCondition(ts, u);CHKERRQ(ierr);
385*c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
386*c4762a1bSJed Brown   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "kinematics", &u);CHKERRQ(ierr);
387*c4762a1bSJed Brown 
388*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
389*c4762a1bSJed Brown   ierr = DMDestroy(&sw);CHKERRQ(ierr);
390*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
391*c4762a1bSJed Brown   ierr = PetscFinalize();
392*c4762a1bSJed Brown   return ierr;
393*c4762a1bSJed Brown }
394*c4762a1bSJed Brown 
395*c4762a1bSJed Brown /*TEST
396*c4762a1bSJed Brown 
397*c4762a1bSJed Brown    build:
398*c4762a1bSJed Brown      requires: triangle !single !complex
399*c4762a1bSJed Brown    test:
400*c4762a1bSJed Brown      suffix: 1
401*c4762a1bSJed Brown      args: -dm_plex_box_faces 1,1 -ts_basicsymplectic_type 1 -ts_convergence_estimate -convest_num_refine 2 -dm_view -sw_view -monitor -output_step 50 -error
402*c4762a1bSJed Brown    test:
403*c4762a1bSJed Brown      suffix: 2
404*c4762a1bSJed Brown      args: -dm_plex_box_faces 1,1 -ts_basicsymplectic_type 2 -ts_convergence_estimate -convest_num_refine 2 -dm_view -sw_view -monitor -output_step 50 -error
405*c4762a1bSJed Brown    test:
406*c4762a1bSJed Brown      suffix: 3
407*c4762a1bSJed Brown      args: -dm_plex_box_faces 1,1 -ts_basicsymplectic_type 3 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -sw_view -monitor -output_step 50 -error
408*c4762a1bSJed Brown    test:
409*c4762a1bSJed Brown      suffix: 4
410*c4762a1bSJed Brown      args: -dm_plex_box_faces 1,1 -ts_basicsymplectic_type 4 -ts_convergence_estimate -convest_num_refine 2 -ts_dt 0.0001 -dm_view -sw_view -monitor -output_step 50 -error
411*c4762a1bSJed Brown 
412*c4762a1bSJed Brown TEST*/
413