xref: /petsc/src/dm/impls/swarm/tests/ex4.c (revision 3157dc651feb7ccf0dc010bd0b250a5e125c1a17)
1*3157dc65SMatthew G. Knepley static char help[] = "Testing integrators on the simple harmonic oscillator\n";
2c4762a1bSJed Brown 
3*3157dc65SMatthew G. Knepley /*
4*3157dc65SMatthew G. Knepley   In order to check long time behavior, you can give
5*3157dc65SMatthew G. Knepley 
6*3157dc65SMatthew G. Knepley      -ts_max_steps 10000 -ts_max_time 10.0 -output_step 1000
7*3157dc65SMatthew G. Knepley 
8*3157dc65SMatthew G. Knepley   To look at the long time behavior for different integrators, we can use the nergy monitor. Below we see that Euler is almost 100 times worse at conserving energy than the symplectic integrator of the same order.
9*3157dc65SMatthew G. Knepley 
10*3157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_1"
11*3157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 10000 -ts_max_time 10.0 -output_step 1000 -ts_type euler" | grep error
12*3157dc65SMatthew G. Knepley 
13*3157dc65SMatthew G. Knepley       energy/exact energy 398.236 / 396.608 (0.4104399231)
14*3157dc65SMatthew G. Knepley       energy/exact energy 1579.52 / 1573.06 (0.4104399231)
15*3157dc65SMatthew G. Knepley       energy/exact energy 397.421 / 396.608 (0.2050098479)
16*3157dc65SMatthew G. Knepley       energy/exact energy 1576.29 / 1573.06 (0.2050098479)
17*3157dc65SMatthew G. Knepley       energy/exact energy 397.014 / 396.608 (0.1024524454)
18*3157dc65SMatthew G. Knepley       energy/exact energy 1574.68 / 1573.06 (0.1024524454)
19*3157dc65SMatthew G. Knepley 
20*3157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_1"
21*3157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 10000 -ts_max_time 10.0 -output_step 1000" | grep error
22*3157dc65SMatthew G. Knepley 
23*3157dc65SMatthew G. Knepley       energy/exact energy 396.579 / 396.608 (0.0074080434)
24*3157dc65SMatthew G. Knepley       energy/exact energy 1572.95 / 1573.06 (0.0074080434)
25*3157dc65SMatthew G. Knepley       energy/exact energy 396.593 / 396.608 (0.0037040885)
26*3157dc65SMatthew G. Knepley       energy/exact energy 1573.01 / 1573.06 (0.0037040886)
27*3157dc65SMatthew G. Knepley       energy/exact energy 396.601 / 396.608 (0.0018520613)
28*3157dc65SMatthew G. Knepley       energy/exact energy 1573.04 / 1573.06 (0.0018520613)
29*3157dc65SMatthew G. Knepley 
30*3157dc65SMatthew G. Knepley   We can look at third order integrators in the same way, but we need to use more steps.
31*3157dc65SMatthew G. Knepley 
32*3157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_1"
33*3157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 1000000 -ts_max_time 1000.0 -output_step 100000 -ts_type rk -ts_adapt_type none" | grep error
34*3157dc65SMatthew G. Knepley 
35*3157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000013981)
36*3157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000013981)
37*3157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000001747)
38*3157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000001748)
39*3157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000218)
40*3157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000218)
41*3157dc65SMatthew G. Knepley 
42*3157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_3"
43*3157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 1000000 -ts_max_time 1000.0 -output_step 100000 -ts_adapt_type none" | grep error
44*3157dc65SMatthew G. Knepley 
45*3157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000007)
46*3157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000007)
47*3157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000001)
48*3157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000001)
49*3157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000000)
50*3157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000000)
51*3157dc65SMatthew G. Knepley */
52*3157dc65SMatthew G. Knepley 
53c4762a1bSJed Brown #include <petscts.h>
54*3157dc65SMatthew G. Knepley #include <petscdmplex.h>
55*3157dc65SMatthew G. Knepley #include <petscdmswarm.h>
56*3157dc65SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>  /* For norm and dot */
57c4762a1bSJed Brown 
58c4762a1bSJed Brown typedef struct {
59c4762a1bSJed Brown   PetscReal omega; /* Oscillation frequency omega */
60c4762a1bSJed Brown   PetscBool error; /* Flag for printing the error */
61c4762a1bSJed Brown   PetscInt  ostep; /* print the energy at each ostep time steps */
62c4762a1bSJed Brown } AppCtx;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   PetscErrorCode ierr;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   PetscFunctionBeginUser;
69c4762a1bSJed Brown   options->omega = 64.0;
70*3157dc65SMatthew G. Knepley   options->error = PETSC_FALSE;
71c4762a1bSJed Brown   options->ostep = 100;
72c4762a1bSJed Brown 
73*3157dc65SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Harmonic Oscillator Options", "DMSWARM");PetscCall(ierr);
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-omega", "Oscillator frequency", "ex4.c", options->omega, &options->omega, PETSC_NULL));
75*3157dc65SMatthew G. Knepley   PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL));
76*3157dc65SMatthew G. Knepley   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, PETSC_NULL));
779566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
78c4762a1bSJed Brown   PetscFunctionReturn(0);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81*3157dc65SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
859566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
869566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
879566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
88c4762a1bSJed Brown   PetscFunctionReturn(0);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91*3157dc65SMatthew G. Knepley static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
92c4762a1bSJed Brown {
93*3157dc65SMatthew G. Knepley   PetscReal v0[1] = {0.}; // Initialize velocities to 0
94*3157dc65SMatthew G. Knepley   PetscInt  dim;
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   PetscFunctionBeginUser;
979566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
989566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), sw));
999566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
1009566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
1019566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1029566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
103*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
104*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
105*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
106*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL));
1079566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
108*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
109*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmInitializeCoordinates(*sw));
110*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0));
1119566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
112*3157dc65SMatthew G. Knepley   PetscCall(DMSetApplicationContext(*sw, user));
1139566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *sw, "Particles"));
1149566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
115*3157dc65SMatthew G. Knepley   {
116*3157dc65SMatthew G. Knepley     Vec gc, gc0;
117*3157dc65SMatthew G. Knepley 
118*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
119*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0));
120*3157dc65SMatthew G. Knepley     PetscCall(VecCopy(gc, gc0));
121*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
122*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0));
123*3157dc65SMatthew G. Knepley   }
124c4762a1bSJed Brown   PetscFunctionReturn(0);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127*3157dc65SMatthew G. Knepley static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
128c4762a1bSJed Brown {
129*3157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *) ctx)->omega;
130*3157dc65SMatthew G. Knepley   DM                 sw;
131c4762a1bSJed Brown   const PetscScalar *u;
132*3157dc65SMatthew G. Knepley   PetscScalar       *g;
133*3157dc65SMatthew G. Knepley   PetscInt           dim, d, Np, p;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   PetscFunctionBeginUser;
136*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
137*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
139*3157dc65SMatthew G. Knepley   PetscCall(VecGetArray(G, &g));
140*3157dc65SMatthew G. Knepley   PetscCall(VecGetArrayRead(U, &u));
141*3157dc65SMatthew G. Knepley   Np /= 2*dim;
142c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
143*3157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
144*3157dc65SMatthew G. Knepley       g[(p*2+0)*dim + d] = u[(p*2+1)*dim + d];
145*3157dc65SMatthew G. Knepley       g[(p*2+1)*dim + d] = -PetscSqr(omega)*u[(p*2+0)*dim+d];
146*3157dc65SMatthew G. Knepley     }
147c4762a1bSJed Brown   }
1489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
149*3157dc65SMatthew G. Knepley   PetscCall(VecRestoreArray(G, &g));
150c4762a1bSJed Brown   PetscFunctionReturn(0);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153*3157dc65SMatthew G. Knepley /* J_{ij} = dF_i/dx_j
154*3157dc65SMatthew G. Knepley    J_p = (  0   1)
155*3157dc65SMatthew G. Knepley          (-w^2  0)
156*3157dc65SMatthew G. Knepley */
157*3157dc65SMatthew G. Knepley static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U , Mat J, Mat P, void *ctx)
158c4762a1bSJed Brown {
159*3157dc65SMatthew G. Knepley   PetscScalar vals[4] = {0., 1., -PetscSqr(((AppCtx *) ctx)->omega), 0.};
160*3157dc65SMatthew G. Knepley   DM          sw;
161*3157dc65SMatthew G. Knepley   PetscInt    dim, d, Np, p, rStart;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   PetscFunctionBeginUser;
164*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
165*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
1669566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
167*3157dc65SMatthew G. Knepley   PetscCall(MatGetOwnershipRange(J, &rStart, NULL));
168*3157dc65SMatthew G. Knepley   Np /= 2*dim;
169c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
170*3157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
171*3157dc65SMatthew G. Knepley       const PetscInt rows[2] = {(p*2+0)*dim+d + rStart, (p*2+1)*dim+d + rStart};
172*3157dc65SMatthew G. Knepley       PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
173c4762a1bSJed Brown     }
174*3157dc65SMatthew G. Knepley   }
175*3157dc65SMatthew G. Knepley   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
176*3157dc65SMatthew G. Knepley   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180*3157dc65SMatthew G. Knepley static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
18140be0ff1SMatthew G. Knepley {
18240be0ff1SMatthew G. Knepley   const PetscScalar *v;
18340be0ff1SMatthew G. Knepley   PetscScalar       *xres;
18440be0ff1SMatthew G. Knepley   PetscInt           Np, p;
18540be0ff1SMatthew G. Knepley 
18640be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Xres, &Np));
188*3157dc65SMatthew G. Knepley   PetscCall(VecGetArrayRead(V, &v));
189*3157dc65SMatthew G. Knepley   PetscCall(VecGetArray(Xres, &xres));
190*3157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) xres[p] = v[p];
1919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
1929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xres, &xres));
19340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
19440be0ff1SMatthew G. Knepley }
19540be0ff1SMatthew G. Knepley 
196*3157dc65SMatthew G. Knepley static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
19740be0ff1SMatthew G. Knepley {
198*3157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *) ctx)->omega;
19940be0ff1SMatthew G. Knepley   const PetscScalar *x;
20040be0ff1SMatthew G. Knepley   PetscScalar       *vres;
20140be0ff1SMatthew G. Knepley   PetscInt           Np, p;
20240be0ff1SMatthew G. Knepley 
20340be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
2049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vres, &vres));
2059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2069566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Vres, &Np));
207*3157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) vres[p] = -PetscSqr(omega)*x[p];
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vres, &vres));
21040be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
21140be0ff1SMatthew G. Knepley }
21240be0ff1SMatthew G. Knepley 
213*3157dc65SMatthew G. Knepley PetscErrorCode RHSJacobianS(TS ts, PetscReal t, Vec U, Mat S, void *ctx)
21440be0ff1SMatthew G. Knepley {
215*3157dc65SMatthew G. Knepley   PetscScalar vals[4] = {0., 1., -1., 0.};
216*3157dc65SMatthew G. Knepley   DM          sw;
217*3157dc65SMatthew G. Knepley   PetscInt    dim, d, Np, p, rStart;
21840be0ff1SMatthew G. Knepley 
21940be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
220*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
221*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
222*3157dc65SMatthew G. Knepley   PetscCall(VecGetLocalSize(U, &Np));
223*3157dc65SMatthew G. Knepley   PetscCall(MatGetOwnershipRange(S, &rStart, NULL));
224*3157dc65SMatthew G. Knepley   Np /= 2*dim;
225*3157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
226*3157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
227*3157dc65SMatthew G. Knepley       const PetscInt rows[2] = {(p*2+0)*dim+d + rStart, (p*2+1)*dim+d + rStart};
228*3157dc65SMatthew G. Knepley       PetscCall(MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES));
229*3157dc65SMatthew G. Knepley     }
230*3157dc65SMatthew G. Knepley   }
231*3157dc65SMatthew G. Knepley   PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY));
232*3157dc65SMatthew G. Knepley   PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY));
233*3157dc65SMatthew G. Knepley   PetscFunctionReturn(0);
234*3157dc65SMatthew G. Knepley }
235*3157dc65SMatthew G. Knepley 
236*3157dc65SMatthew G. Knepley PetscErrorCode RHSObjectiveF(TS ts, PetscReal t, Vec U, PetscScalar *F, void *ctx)
237*3157dc65SMatthew G. Knepley {
238*3157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *) ctx)->omega;
239*3157dc65SMatthew G. Knepley   DM                 sw;
240*3157dc65SMatthew G. Knepley   const PetscScalar *u;
241*3157dc65SMatthew G. Knepley   PetscInt           dim, Np, p;
242*3157dc65SMatthew G. Knepley 
243*3157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
244*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
245*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
246*3157dc65SMatthew G. Knepley   PetscCall(VecGetArrayRead(U, &u));
247*3157dc65SMatthew G. Knepley   PetscCall(VecGetLocalSize(U, &Np));
248*3157dc65SMatthew G. Knepley   Np /= 2*dim;
249*3157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
250*3157dc65SMatthew G. Knepley     const PetscReal x2 = DMPlex_DotRealD_Internal(dim, &u[(p*2+0)*dim], &u[(p*2+0)*dim]);
251*3157dc65SMatthew G. Knepley     const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p*2+1)*dim], &u[(p*2+1)*dim]);
252*3157dc65SMatthew G. Knepley     const PetscReal E  = 0.5*(v2 + PetscSqr(omega)*x2);
253*3157dc65SMatthew G. Knepley 
254*3157dc65SMatthew G. Knepley     *F += E;
255*3157dc65SMatthew G. Knepley   }
256*3157dc65SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(U, &u));
257*3157dc65SMatthew G. Knepley   PetscFunctionReturn(0);
258*3157dc65SMatthew G. Knepley }
259*3157dc65SMatthew G. Knepley 
260*3157dc65SMatthew G. Knepley /* dF/dx = omega^2 x   dF/dv = v */
261*3157dc65SMatthew G. Knepley PetscErrorCode RHSFunctionG(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
262*3157dc65SMatthew G. Knepley {
263*3157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *) ctx)->omega;
264*3157dc65SMatthew G. Knepley   DM                 sw;
265*3157dc65SMatthew G. Knepley   const PetscScalar *u;
266*3157dc65SMatthew G. Knepley   PetscScalar       *g;
267*3157dc65SMatthew G. Knepley   PetscInt           dim, d, Np, p;
268*3157dc65SMatthew G. Knepley 
269*3157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
270*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
271*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
2739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2749566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
275*3157dc65SMatthew G. Knepley   Np /= 2*dim;
27640be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
277*3157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
278*3157dc65SMatthew G. Knepley       g[(p*2+0)*dim + d] = PetscSqr(omega)*u[(p*2+0)*dim+d];
279*3157dc65SMatthew G. Knepley       g[(p*2+1)*dim + d] = u[(p*2+1)*dim + d];
280*3157dc65SMatthew G. Knepley     }
28140be0ff1SMatthew G. Knepley   }
2829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
28440be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
28540be0ff1SMatthew G. Knepley }
28640be0ff1SMatthew G. Knepley 
287*3157dc65SMatthew G. Knepley static PetscErrorCode CreateSolution(TS ts)
28840be0ff1SMatthew G. Knepley {
289*3157dc65SMatthew G. Knepley   DM       sw;
290*3157dc65SMatthew G. Knepley   Vec      u;
291*3157dc65SMatthew G. Knepley   PetscInt dim, Np;
29240be0ff1SMatthew G. Knepley 
293*3157dc65SMatthew G. Knepley   PetscFunctionBegin;
294*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
295*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
296*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(sw, &Np));
297*3157dc65SMatthew G. Knepley   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
298*3157dc65SMatthew G. Knepley   PetscCall(VecSetBlockSize(u, dim));
299*3157dc65SMatthew G. Knepley   PetscCall(VecSetSizes(u, 2*Np*dim, PETSC_DECIDE));
300*3157dc65SMatthew G. Knepley   PetscCall(VecSetUp(u));
301*3157dc65SMatthew G. Knepley   PetscCall(TSSetSolution(ts, u));
302*3157dc65SMatthew G. Knepley   PetscCall(VecDestroy(&u));
30340be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
30440be0ff1SMatthew G. Knepley }
305db4653c2SDaniel Finn 
306*3157dc65SMatthew G. Knepley static PetscErrorCode SetProblem(TS ts)
30740be0ff1SMatthew G. Knepley {
308*3157dc65SMatthew G. Knepley   AppCtx *user;
309*3157dc65SMatthew G. Knepley   DM      sw;
31040be0ff1SMatthew G. Knepley 
311*3157dc65SMatthew G. Knepley   PetscFunctionBegin;
312*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
313*3157dc65SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, (void **) &user));
314*3157dc65SMatthew G. Knepley   // Define unified system for (X, V)
315*3157dc65SMatthew G. Knepley   {
316*3157dc65SMatthew G. Knepley     Mat      J;
317*3157dc65SMatthew G. Knepley     PetscInt dim, Np;
318*3157dc65SMatthew G. Knepley 
319*3157dc65SMatthew G. Knepley     PetscCall(DMGetDimension(sw, &dim));
320*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &Np));
321*3157dc65SMatthew G. Knepley     PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
322*3157dc65SMatthew G. Knepley     PetscCall(MatSetSizes(J, 2*Np*dim, 2*Np*dim, PETSC_DECIDE, PETSC_DECIDE));
323*3157dc65SMatthew G. Knepley     PetscCall(MatSetBlockSize(J, 2*dim));
324*3157dc65SMatthew G. Knepley     PetscCall(MatSetFromOptions(J));
325*3157dc65SMatthew G. Knepley     PetscCall(MatSetUp(J));
326*3157dc65SMatthew G. Knepley     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user));
327*3157dc65SMatthew G. Knepley     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user));
328*3157dc65SMatthew G. Knepley     PetscCall(MatDestroy(&J));
32940be0ff1SMatthew G. Knepley   }
330*3157dc65SMatthew G. Knepley   // Define split system for X and V
331*3157dc65SMatthew G. Knepley   {
332*3157dc65SMatthew G. Knepley     IS              isx, isv, istmp;
333*3157dc65SMatthew G. Knepley     const PetscInt *idx;
334*3157dc65SMatthew G. Knepley     PetscInt        dim, Np;
335*3157dc65SMatthew G. Knepley 
336*3157dc65SMatthew G. Knepley     PetscCall(DMGetDimension(sw, &dim));
337*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &Np));
338*3157dc65SMatthew G. Knepley     PetscCall(ISCreateStride(PETSC_COMM_SELF, Np, 0, 2, &istmp));
339*3157dc65SMatthew G. Knepley     PetscCall(ISGetIndices(istmp, &idx));
340*3157dc65SMatthew G. Knepley     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, Np, idx, PETSC_COPY_VALUES, &isx));
341*3157dc65SMatthew G. Knepley     PetscCall(ISRestoreIndices(istmp, &idx));
342*3157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&istmp));
343*3157dc65SMatthew G. Knepley     PetscCall(ISCreateStride(PETSC_COMM_SELF, Np, 1, 2, &istmp));
344*3157dc65SMatthew G. Knepley     PetscCall(ISGetIndices(istmp, &idx));
345*3157dc65SMatthew G. Knepley     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, Np, idx, PETSC_COPY_VALUES, &isv));
346*3157dc65SMatthew G. Knepley     PetscCall(ISRestoreIndices(istmp, &idx));
347*3157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&istmp));
348*3157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetIS(ts, "position", isx));
349*3157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetIS(ts, "momentum", isv));
350*3157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&isx));
351*3157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&isv));
352*3157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user));
353*3157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user));
354*3157dc65SMatthew G. Knepley   }
355*3157dc65SMatthew G. Knepley   // Define symplectic formulation U_t = S . G, where G = grad F
356*3157dc65SMatthew G. Knepley   {
357*3157dc65SMatthew G. Knepley     PetscCall(TSDiscGradSetFormulation(ts, RHSJacobianS, RHSObjectiveF, RHSFunctionG, user));
358*3157dc65SMatthew G. Knepley   }
35940be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
36040be0ff1SMatthew G. Knepley }
36140be0ff1SMatthew G. Knepley 
362*3157dc65SMatthew G. Knepley static PetscErrorCode InitializeSolve(TS ts, Vec u)
36340be0ff1SMatthew G. Knepley {
364*3157dc65SMatthew G. Knepley   DM      sw;
365*3157dc65SMatthew G. Knepley   Vec     gc, gc0;
366*3157dc65SMatthew G. Knepley   IS      isx, isv;
367*3157dc65SMatthew G. Knepley   AppCtx *user;
36840be0ff1SMatthew G. Knepley 
36940be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
370*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
371*3157dc65SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, &user));
372*3157dc65SMatthew G. Knepley   {
373*3157dc65SMatthew G. Knepley     PetscReal v0[1] = {1.};
374db4653c2SDaniel Finn 
375*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmInitializeCoordinates(sw));
376*3157dc65SMatthew G. Knepley     PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0));
377*3157dc65SMatthew G. Knepley     PetscCall(SetProblem(ts));
378*3157dc65SMatthew G. Knepley   }
379*3157dc65SMatthew G. Knepley   PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
380*3157dc65SMatthew G. Knepley   PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
381*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
382*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates",    &gc0));
383*3157dc65SMatthew G. Knepley   PetscCall(VecCopy(gc, gc0));
384*3157dc65SMatthew G. Knepley   PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc));
385*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
386*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates",    &gc0));
387*3157dc65SMatthew G. Knepley   PetscCall(VecISSet(u, isv, 0.));
388*3157dc65SMatthew G. Knepley   PetscFunctionReturn(0);
389*3157dc65SMatthew G. Knepley }
390*3157dc65SMatthew G. Knepley 
391*3157dc65SMatthew G. Knepley static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
392*3157dc65SMatthew G. Knepley {
393*3157dc65SMatthew G. Knepley   MPI_Comm           comm;
394*3157dc65SMatthew G. Knepley   DM                 sw;
395*3157dc65SMatthew G. Knepley   AppCtx            *user;
396*3157dc65SMatthew G. Knepley   const PetscScalar *u;
397*3157dc65SMatthew G. Knepley   const PetscReal   *coords;
398*3157dc65SMatthew G. Knepley   PetscScalar       *e;
399*3157dc65SMatthew G. Knepley   PetscReal          t;
400*3157dc65SMatthew G. Knepley   PetscInt           dim, d, Np, p;
401*3157dc65SMatthew G. Knepley 
402*3157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
403*3157dc65SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject) ts, &comm));
404*3157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
405*3157dc65SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, &user));
406*3157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
407*3157dc65SMatthew G. Knepley   PetscCall(TSGetSolveTime(ts, &t));
408*3157dc65SMatthew G. Knepley   PetscCall(VecGetArray(E, &e));
4099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
4109566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
411*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **) &coords));
412*3157dc65SMatthew G. Knepley   Np /= 2*dim;
41340be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
414*3157dc65SMatthew G. Knepley     const PetscReal omega = user->omega;
415*3157dc65SMatthew G. Knepley     const PetscReal ct = PetscCosReal(omega*t);
416*3157dc65SMatthew G. Knepley     const PetscReal st = PetscSinReal(omega*t);
417*3157dc65SMatthew G. Knepley     const PetscReal x0 = DMPlex_NormD_Internal(dim, &coords[p*dim]);
418*3157dc65SMatthew G. Knepley     const PetscReal ex =  x0*ct;
419*3157dc65SMatthew G. Knepley     const PetscReal ev = -x0*omega*st;
420*3157dc65SMatthew G. Knepley     const PetscReal x  = DMPlex_DotRealD_Internal(dim, &u[(p*2+0)*dim], &coords[p*dim]) / x0;
421*3157dc65SMatthew G. Knepley     const PetscReal v  = DMPlex_DotRealD_Internal(dim, &u[(p*2+1)*dim], &coords[p*dim]) / x0;
422*3157dc65SMatthew G. Knepley 
423*3157dc65SMatthew G. Knepley     if (user->error) {
424*3157dc65SMatthew G. Knepley       const PetscReal en   = 0.5*(v*v + PetscSqr(omega)*x*x);
425*3157dc65SMatthew G. Knepley       const PetscReal exen = 0.5*PetscSqr(omega*x0);
426*3157dc65SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "p%D error [%.2g %.2g] sol [%.6lf %.6lf] exact [%.6lf %.6lf] energy/exact energy %g / %g (%.10lf%%)\n", p, (double) PetscAbsReal(x-ex), (double) PetscAbsReal(v-ev), (double) x, (double) v, (double) ex, (double) ev, (double) en, (double) exen, (double) PetscAbsReal(exen - en)*100./exen));
427*3157dc65SMatthew G. Knepley     }
428*3157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
429*3157dc65SMatthew G. Knepley       e[(p*2+0)*dim+d] = u[(p*2+0)*dim+d] - coords[p*dim+d]*ct;
430*3157dc65SMatthew G. Knepley       e[(p*2+1)*dim+d] = u[(p*2+1)*dim+d] + coords[p*dim+d]*omega*st;
431*3157dc65SMatthew G. Knepley     }
432*3157dc65SMatthew G. Knepley   }
433*3157dc65SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **) &coords));
434*3157dc65SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(U, &u));
435*3157dc65SMatthew G. Knepley   PetscCall(VecRestoreArray(E, &e));
436*3157dc65SMatthew G. Knepley   PetscFunctionReturn(0);
437*3157dc65SMatthew G. Knepley }
438*3157dc65SMatthew G. Knepley 
439*3157dc65SMatthew G. Knepley static PetscErrorCode EnergyMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
440*3157dc65SMatthew G. Knepley {
441*3157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *) ctx)->omega;
442*3157dc65SMatthew G. Knepley   const PetscInt     ostep = ((AppCtx *) ctx)->ostep;
443*3157dc65SMatthew G. Knepley   DM                 sw;
444*3157dc65SMatthew G. Knepley   const PetscScalar *u;
445*3157dc65SMatthew G. Knepley   PetscReal          dt;
446*3157dc65SMatthew G. Knepley   PetscInt           dim, Np, p;
447*3157dc65SMatthew G. Knepley   MPI_Comm           comm;
448*3157dc65SMatthew G. Knepley 
449*3157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
450*3157dc65SMatthew G. Knepley   if (step%ostep == 0) {
451*3157dc65SMatthew G. Knepley     PetscCall(PetscObjectGetComm((PetscObject) ts, &comm));
452*3157dc65SMatthew G. Knepley     PetscCall(TSGetDM(ts, &sw));
453*3157dc65SMatthew G. Knepley     PetscCall(TSGetTimeStep(ts, &dt));
454*3157dc65SMatthew G. Knepley     PetscCall(DMGetDimension(sw, &dim));
455*3157dc65SMatthew G. Knepley     PetscCall(VecGetArrayRead(U, &u));
456*3157dc65SMatthew G. Knepley     PetscCall(VecGetLocalSize(U, &Np));
457*3157dc65SMatthew G. Knepley     Np /= 2*dim;
458*3157dc65SMatthew G. Knepley     if (!step) PetscCall(PetscPrintf(comm, "Time     Step Part     Energy Mod Energy\n"));
459*3157dc65SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
460*3157dc65SMatthew G. Knepley       const PetscReal x2 = DMPlex_DotRealD_Internal(dim, &u[(p*2+0)*dim], &u[(p*2+0)*dim]);
461*3157dc65SMatthew G. Knepley       const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p*2+1)*dim], &u[(p*2+1)*dim]);
462*3157dc65SMatthew G. Knepley       const PetscReal E  = 0.5*(v2 + PetscSqr(omega)*x2);
463*3157dc65SMatthew G. Knepley       const PetscReal mE = 0.5*(v2 + PetscSqr(omega)*x2 - PetscSqr(omega)*dt*PetscSqrtReal(x2*v2));
464*3157dc65SMatthew G. Knepley 
465*3157dc65SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "%.6lf %4D %4D %10.4lf %10.4lf\n", t, step, p, (double) E, (double) mE));
46640be0ff1SMatthew G. Knepley     }
4679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U, &u));
468*3157dc65SMatthew G. Knepley   }
46940be0ff1SMatthew G. Knepley   PetscFunctionReturn(0);
47040be0ff1SMatthew G. Knepley }
47140be0ff1SMatthew G. Knepley 
472c4762a1bSJed Brown int main(int argc, char **argv)
473c4762a1bSJed Brown {
474*3157dc65SMatthew G. Knepley   DM     dm, sw;
475*3157dc65SMatthew G. Knepley   TS     ts;
476*3157dc65SMatthew G. Knepley   Vec    u;
477c4762a1bSJed Brown   AppCtx user;
478c4762a1bSJed Brown 
4799566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
480*3157dc65SMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
481*3157dc65SMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
482*3157dc65SMatthew G. Knepley   PetscCall(CreateSwarm(dm, &user, &sw));
4839566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(sw, &user));
484c4762a1bSJed Brown 
485*3157dc65SMatthew G. Knepley   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
4869566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
4879566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, sw));
4889566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 0.1));
4899566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.00001));
4909566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 100));
4919566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
492*3157dc65SMatthew G. Knepley   PetscCall(TSMonitorSet(ts, EnergyMonitor, &user, NULL));
4939566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
4949566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
4959566063dSJacob Faibussowitsch   PetscCall(TSSetComputeExactError(ts, ComputeError));
49640be0ff1SMatthew G. Knepley 
497*3157dc65SMatthew G. Knepley   PetscCall(CreateSolution(ts));
498*3157dc65SMatthew G. Knepley   PetscCall(TSGetSolution(ts, &u));
4999566063dSJacob Faibussowitsch   PetscCall(TSComputeInitialCondition(ts, u));
500*3157dc65SMatthew G. Knepley   PetscCall(TSSolve(ts, NULL));
501c4762a1bSJed Brown 
5029566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
5039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
5049566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
5059566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
506b122ec5aSJacob Faibussowitsch   return 0;
507c4762a1bSJed Brown }
508c4762a1bSJed Brown 
509c4762a1bSJed Brown /*TEST
510c4762a1bSJed Brown 
511c4762a1bSJed Brown    build:
512c4762a1bSJed Brown      requires: triangle !single !complex
51340be0ff1SMatthew G. Knepley 
514*3157dc65SMatthew G. Knepley    testset:
515*3157dc65SMatthew G. Knepley      args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1 \
516*3157dc65SMatthew G. Knepley            -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
517*3157dc65SMatthew G. Knepley            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
518*3157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
519c4762a1bSJed Brown      test:
520*3157dc65SMatthew G. Knepley        suffix: bsi_1
521*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 1
522*3157dc65SMatthew G. Knepley      test:
523*3157dc65SMatthew G. Knepley        suffix: bsi_2
524*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 2
525*3157dc65SMatthew G. Knepley      test:
526*3157dc65SMatthew G. Knepley        suffix: bsi_3
527*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 3
528*3157dc65SMatthew G. Knepley      test:
529*3157dc65SMatthew G. Knepley        suffix: bsi_4
530*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
53140be0ff1SMatthew G. Knepley 
532*3157dc65SMatthew G. Knepley    testset:
533*3157dc65SMatthew G. Knepley      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1 \
534*3157dc65SMatthew G. Knepley            -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
535*3157dc65SMatthew G. Knepley            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
536*3157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
537c4762a1bSJed Brown      test:
538*3157dc65SMatthew G. Knepley        suffix: bsi_2d_1
539*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 1
540*3157dc65SMatthew G. Knepley      test:
541*3157dc65SMatthew G. Knepley        suffix: bsi_2d_2
542*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 2
543*3157dc65SMatthew G. Knepley      test:
544*3157dc65SMatthew G. Knepley        suffix: bsi_2d_3
545*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 3
546*3157dc65SMatthew G. Knepley      test:
547*3157dc65SMatthew G. Knepley        suffix: bsi_2d_4
548*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
54940be0ff1SMatthew G. Knepley 
550*3157dc65SMatthew G. Knepley    testset:
551*3157dc65SMatthew G. Knepley      args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1,-1 -dm_plex_box_upper 1,1,1 \
552*3157dc65SMatthew G. Knepley            -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
553*3157dc65SMatthew G. Knepley            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
554*3157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
555c4762a1bSJed Brown      test:
556*3157dc65SMatthew G. Knepley        suffix: bsi_3d_1
557*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 1
558*3157dc65SMatthew G. Knepley      test:
559*3157dc65SMatthew G. Knepley        suffix: bsi_3d_2
560*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 2
561*3157dc65SMatthew G. Knepley      test:
562*3157dc65SMatthew G. Knepley        suffix: bsi_3d_3
563*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 3
564*3157dc65SMatthew G. Knepley      test:
565*3157dc65SMatthew G. Knepley        suffix: bsi_3d_4
566*3157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
56740be0ff1SMatthew G. Knepley 
568*3157dc65SMatthew G. Knepley    testset:
569*3157dc65SMatthew G. Knepley      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
570*3157dc65SMatthew G. Knepley            -ts_type theta -ts_theta_theta 0.5 -ts_convergence_estimate -convest_num_refine 2 \
571*3157dc65SMatthew G. Knepley              -mat_type baij -ksp_error_if_not_converged -pc_type lu \
572*3157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
57340be0ff1SMatthew G. Knepley      test:
574*3157dc65SMatthew G. Knepley        suffix: im_1d
575*3157dc65SMatthew G. Knepley        args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1
576*3157dc65SMatthew G. Knepley      test:
577*3157dc65SMatthew G. Knepley        suffix: im_2d
578*3157dc65SMatthew G. Knepley        args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1
579*3157dc65SMatthew G. Knepley      test:
580*3157dc65SMatthew G. Knepley        suffix: im_3d
581*3157dc65SMatthew G. Knepley        args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1,-1,-1 -dm_plex_box_upper 1,1,1
58240be0ff1SMatthew G. Knepley 
583*3157dc65SMatthew G. Knepley    testset:
584*3157dc65SMatthew G. Knepley      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
585*3157dc65SMatthew G. Knepley            -ts_type discgrad -ts_discgrad_gonzalez -ts_convergence_estimate -convest_num_refine 2 \
586*3157dc65SMatthew G. Knepley              -mat_type baij -ksp_error_if_not_converged -pc_type lu \
587*3157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
58840be0ff1SMatthew G. Knepley      test:
589*3157dc65SMatthew G. Knepley        suffix: dg_1d
590*3157dc65SMatthew G. Knepley        args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1
591db4653c2SDaniel Finn      test:
592*3157dc65SMatthew G. Knepley        suffix: dg_2d
593*3157dc65SMatthew G. Knepley        args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -1,-1 -dm_plex_box_upper 1,1
594*3157dc65SMatthew G. Knepley      test:
595*3157dc65SMatthew G. Knepley        suffix: dg_3d
596*3157dc65SMatthew G. Knepley        args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_faces 1,1,1 -dm_plex_box_lower -1,-1,-1 -dm_plex_box_upper 1,1,1
597c4762a1bSJed Brown 
598c4762a1bSJed Brown TEST*/
599