xref: /petsc/src/dm/impls/swarm/tests/ex4.c (revision f940b0e3a9f0ece732f2385619a7ea1fae51a73a)
13157dc65SMatthew G. Knepley static char help[] = "Testing integrators on the simple harmonic oscillator\n";
2c4762a1bSJed Brown 
33157dc65SMatthew G. Knepley /*
43157dc65SMatthew G. Knepley   In order to check long time behavior, you can give
53157dc65SMatthew G. Knepley 
63157dc65SMatthew G. Knepley      -ts_max_steps 10000 -ts_max_time 10.0 -output_step 1000
73157dc65SMatthew G. Knepley 
83157dc65SMatthew 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.
93157dc65SMatthew G. Knepley 
103157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_1"
113157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 10000 -ts_max_time 10.0 -output_step 1000 -ts_type euler" | grep error
123157dc65SMatthew G. Knepley 
133157dc65SMatthew G. Knepley       energy/exact energy 398.236 / 396.608 (0.4104399231)
143157dc65SMatthew G. Knepley       energy/exact energy 1579.52 / 1573.06 (0.4104399231)
153157dc65SMatthew G. Knepley       energy/exact energy 397.421 / 396.608 (0.2050098479)
163157dc65SMatthew G. Knepley       energy/exact energy 1576.29 / 1573.06 (0.2050098479)
173157dc65SMatthew G. Knepley       energy/exact energy 397.014 / 396.608 (0.1024524454)
183157dc65SMatthew G. Knepley       energy/exact energy 1574.68 / 1573.06 (0.1024524454)
193157dc65SMatthew G. Knepley 
203157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_1"
213157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 10000 -ts_max_time 10.0 -output_step 1000" | grep error
223157dc65SMatthew G. Knepley 
233157dc65SMatthew G. Knepley       energy/exact energy 396.579 / 396.608 (0.0074080434)
243157dc65SMatthew G. Knepley       energy/exact energy 1572.95 / 1573.06 (0.0074080434)
253157dc65SMatthew G. Knepley       energy/exact energy 396.593 / 396.608 (0.0037040885)
263157dc65SMatthew G. Knepley       energy/exact energy 1573.01 / 1573.06 (0.0037040886)
273157dc65SMatthew G. Knepley       energy/exact energy 396.601 / 396.608 (0.0018520613)
283157dc65SMatthew G. Knepley       energy/exact energy 1573.04 / 1573.06 (0.0018520613)
293157dc65SMatthew G. Knepley 
303157dc65SMatthew G. Knepley   We can look at third order integrators in the same way, but we need to use more steps.
313157dc65SMatthew G. Knepley 
323157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_1"
333157dc65SMatthew 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
343157dc65SMatthew G. Knepley 
353157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000013981)
363157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000013981)
373157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000001747)
383157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000001748)
393157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000218)
403157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000218)
413157dc65SMatthew G. Knepley 
423157dc65SMatthew G. Knepley     make -f ./gmakefile test search="dm_impls_swarm_tests-ex4_bsi_3"
433157dc65SMatthew G. Knepley       EXTRA_OPTIONS="-ts_max_steps 1000000 -ts_max_time 1000.0 -output_step 100000 -ts_adapt_type none" | grep error
443157dc65SMatthew G. Knepley 
453157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000007)
463157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000007)
473157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000001)
483157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000001)
493157dc65SMatthew G. Knepley       energy/exact energy 396.608 / 396.608 (0.0000000000)
503157dc65SMatthew G. Knepley       energy/exact energy 1573.06 / 1573.06 (0.0000000000)
513157dc65SMatthew G. Knepley */
523157dc65SMatthew G. Knepley 
53c4762a1bSJed Brown #include <petscts.h>
543157dc65SMatthew G. Knepley #include <petscdmplex.h>
553157dc65SMatthew G. Knepley #include <petscdmswarm.h>
563157dc65SMatthew 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 
64d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
65d71ae5a4SJacob Faibussowitsch {
66c4762a1bSJed Brown   PetscFunctionBeginUser;
67c4762a1bSJed Brown   options->omega = 64.0;
683157dc65SMatthew G. Knepley   options->error = PETSC_FALSE;
69c4762a1bSJed Brown   options->ostep = 100;
70c4762a1bSJed Brown 
71d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Harmonic Oscillator Options", "DMSWARM");
72f3fa974cSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-omega", "Oscillator frequency", "ex4.c", options->omega, &options->omega, NULL));
733157dc65SMatthew G. Knepley   PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex4.c", options->error, &options->error, NULL));
74f3fa974cSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex4.c", options->ostep, &options->ostep, NULL));
75d0609cedSBarry Smith   PetscOptionsEnd();
763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
77c4762a1bSJed Brown }
78c4762a1bSJed Brown 
79d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
80d71ae5a4SJacob Faibussowitsch {
81c4762a1bSJed Brown   PetscFunctionBeginUser;
829566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
839566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
849566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
859566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
90d71ae5a4SJacob Faibussowitsch {
913157dc65SMatthew G. Knepley   PetscReal v0[1] = {0.}; // Initialize velocities to 0
923157dc65SMatthew G. Knepley   PetscInt  dim;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   PetscFunctionBeginUser;
959566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
969566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
979566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
989566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
999566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1009566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
1013157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1023157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
1033157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
1043157dc65SMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL));
1059566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
1063157dc65SMatthew G. Knepley   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
1073157dc65SMatthew G. Knepley   PetscCall(DMSwarmInitializeCoordinates(*sw));
1083157dc65SMatthew G. Knepley   PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0));
1099566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
1103157dc65SMatthew G. Knepley   PetscCall(DMSetApplicationContext(*sw, user));
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
1129566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
1133157dc65SMatthew G. Knepley   {
1143157dc65SMatthew G. Knepley     Vec gc, gc0;
1153157dc65SMatthew G. Knepley 
1163157dc65SMatthew G. Knepley     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
1173157dc65SMatthew G. Knepley     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0));
1183157dc65SMatthew G. Knepley     PetscCall(VecCopy(gc, gc0));
1193157dc65SMatthew G. Knepley     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
1203157dc65SMatthew G. Knepley     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0));
1213157dc65SMatthew G. Knepley   }
122*f940b0e3Sdanofinn   {
123*f940b0e3Sdanofinn     const char *fieldnames[2] = {DMSwarmPICField_coor, "velocity"};
124*f940b0e3Sdanofinn     PetscCall(DMSwarmVectorDefineFields(*sw, 2, fieldnames));
125*f940b0e3Sdanofinn   }
1263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
130d71ae5a4SJacob Faibussowitsch {
1313157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *)ctx)->omega;
1323157dc65SMatthew G. Knepley   DM                 sw;
133c4762a1bSJed Brown   const PetscScalar *u;
1343157dc65SMatthew G. Knepley   PetscScalar       *g;
1353157dc65SMatthew G. Knepley   PetscInt           dim, d, Np, p;
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   PetscFunctionBeginUser;
1383157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
1393157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
1409566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1413157dc65SMatthew G. Knepley   PetscCall(VecGetArray(G, &g));
1423157dc65SMatthew G. Knepley   PetscCall(VecGetArrayRead(U, &u));
1433157dc65SMatthew G. Knepley   Np /= 2 * dim;
144c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
1453157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1463157dc65SMatthew G. Knepley       g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d];
1473157dc65SMatthew G. Knepley       g[(p * 2 + 1) * dim + d] = -PetscSqr(omega) * u[(p * 2 + 0) * dim + d];
1483157dc65SMatthew G. Knepley     }
149c4762a1bSJed Brown   }
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1513157dc65SMatthew G. Knepley   PetscCall(VecRestoreArray(G, &g));
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
1553157dc65SMatthew G. Knepley /* J_{ij} = dF_i/dx_j
1563157dc65SMatthew G. Knepley    J_p = (  0   1)
1573157dc65SMatthew G. Knepley          (-w^2  0)
1583157dc65SMatthew G. Knepley */
159d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, void *ctx)
160d71ae5a4SJacob Faibussowitsch {
1613157dc65SMatthew G. Knepley   PetscScalar vals[4] = {0., 1., -PetscSqr(((AppCtx *)ctx)->omega), 0.};
1623157dc65SMatthew G. Knepley   DM          sw;
1633157dc65SMatthew G. Knepley   PetscInt    dim, d, Np, p, rStart;
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   PetscFunctionBeginUser;
1663157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
1673157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
1689566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
1693157dc65SMatthew G. Knepley   PetscCall(MatGetOwnershipRange(J, &rStart, NULL));
1703157dc65SMatthew G. Knepley   Np /= 2 * dim;
171c4762a1bSJed Brown   for (p = 0; p < Np; ++p) {
1723157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1733157dc65SMatthew G. Knepley       const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
1743157dc65SMatthew G. Knepley       PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
175c4762a1bSJed Brown     }
1763157dc65SMatthew G. Knepley   }
1773157dc65SMatthew G. Knepley   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1783157dc65SMatthew G. Knepley   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
183d71ae5a4SJacob Faibussowitsch {
18440be0ff1SMatthew G. Knepley   const PetscScalar *v;
18540be0ff1SMatthew G. Knepley   PetscScalar       *xres;
18640be0ff1SMatthew G. Knepley   PetscInt           Np, p;
18740be0ff1SMatthew G. Knepley 
18840be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
1899566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Xres, &Np));
1903157dc65SMatthew G. Knepley   PetscCall(VecGetArrayRead(V, &v));
1913157dc65SMatthew G. Knepley   PetscCall(VecGetArray(Xres, &xres));
1923157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) xres[p] = v[p];
1939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
1949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xres, &xres));
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19640be0ff1SMatthew G. Knepley }
19740be0ff1SMatthew G. Knepley 
198d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
199d71ae5a4SJacob Faibussowitsch {
2003157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *)ctx)->omega;
20140be0ff1SMatthew G. Knepley   const PetscScalar *x;
20240be0ff1SMatthew G. Knepley   PetscScalar       *vres;
20340be0ff1SMatthew G. Knepley   PetscInt           Np, p;
20440be0ff1SMatthew G. Knepley 
20540be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vres, &vres));
2079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2089566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(Vres, &Np));
2093157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) vres[p] = -PetscSqr(omega) * x[p];
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vres, &vres));
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21340be0ff1SMatthew G. Knepley }
21440be0ff1SMatthew G. Knepley 
215d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobianS(TS ts, PetscReal t, Vec U, Mat S, void *ctx)
216d71ae5a4SJacob Faibussowitsch {
2173157dc65SMatthew G. Knepley   PetscScalar vals[4] = {0., 1., -1., 0.};
2183157dc65SMatthew G. Knepley   DM          sw;
2193157dc65SMatthew G. Knepley   PetscInt    dim, d, Np, p, rStart;
22040be0ff1SMatthew G. Knepley 
22140be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
2223157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
2233157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2243157dc65SMatthew G. Knepley   PetscCall(VecGetLocalSize(U, &Np));
2253157dc65SMatthew G. Knepley   PetscCall(MatGetOwnershipRange(S, &rStart, NULL));
2263157dc65SMatthew G. Knepley   Np /= 2 * dim;
2273157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2283157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
2293157dc65SMatthew G. Knepley       const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
2303157dc65SMatthew G. Knepley       PetscCall(MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES));
2313157dc65SMatthew G. Knepley     }
2323157dc65SMatthew G. Knepley   }
2333157dc65SMatthew G. Knepley   PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY));
2343157dc65SMatthew G. Knepley   PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY));
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2363157dc65SMatthew G. Knepley }
2373157dc65SMatthew G. Knepley 
238d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSObjectiveF(TS ts, PetscReal t, Vec U, PetscScalar *F, void *ctx)
239d71ae5a4SJacob Faibussowitsch {
2403157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *)ctx)->omega;
2413157dc65SMatthew G. Knepley   DM                 sw;
2423157dc65SMatthew G. Knepley   const PetscScalar *u;
2433157dc65SMatthew G. Knepley   PetscInt           dim, Np, p;
2443157dc65SMatthew G. Knepley 
2453157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
2463157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
2473157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2483157dc65SMatthew G. Knepley   PetscCall(VecGetArrayRead(U, &u));
2493157dc65SMatthew G. Knepley   PetscCall(VecGetLocalSize(U, &Np));
2503157dc65SMatthew G. Knepley   Np /= 2 * dim;
2513157dc65SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2523157dc65SMatthew G. Knepley     const PetscReal x2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 0) * dim], &u[(p * 2 + 0) * dim]);
2533157dc65SMatthew G. Knepley     const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]);
2543157dc65SMatthew G. Knepley     const PetscReal E  = 0.5 * (v2 + PetscSqr(omega) * x2);
2553157dc65SMatthew G. Knepley 
2563157dc65SMatthew G. Knepley     *F += E;
2573157dc65SMatthew G. Knepley   }
2583157dc65SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(U, &u));
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2603157dc65SMatthew G. Knepley }
2613157dc65SMatthew G. Knepley 
2623157dc65SMatthew G. Knepley /* dF/dx = omega^2 x   dF/dv = v */
263d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunctionG(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
264d71ae5a4SJacob Faibussowitsch {
2653157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *)ctx)->omega;
2663157dc65SMatthew G. Knepley   DM                 sw;
2673157dc65SMatthew G. Knepley   const PetscScalar *u;
2683157dc65SMatthew G. Knepley   PetscScalar       *g;
2693157dc65SMatthew G. Knepley   PetscInt           dim, d, Np, p;
2703157dc65SMatthew G. Knepley 
2713157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
2723157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
2733157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
2759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2769566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
2773157dc65SMatthew G. Knepley   Np /= 2 * dim;
27840be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
2793157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
2803157dc65SMatthew G. Knepley       g[(p * 2 + 0) * dim + d] = PetscSqr(omega) * u[(p * 2 + 0) * dim + d];
2813157dc65SMatthew G. Knepley       g[(p * 2 + 1) * dim + d] = u[(p * 2 + 1) * dim + d];
2823157dc65SMatthew G. Knepley     }
28340be0ff1SMatthew G. Knepley   }
2849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28740be0ff1SMatthew G. Knepley }
28840be0ff1SMatthew G. Knepley 
289d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSolution(TS ts)
290d71ae5a4SJacob Faibussowitsch {
2913157dc65SMatthew G. Knepley   DM       sw;
2923157dc65SMatthew G. Knepley   Vec      u;
2933157dc65SMatthew G. Knepley   PetscInt dim, Np;
29440be0ff1SMatthew G. Knepley 
2953157dc65SMatthew G. Knepley   PetscFunctionBegin;
2963157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
2973157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2983157dc65SMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(sw, &Np));
2993157dc65SMatthew G. Knepley   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
3003157dc65SMatthew G. Knepley   PetscCall(VecSetBlockSize(u, dim));
3013157dc65SMatthew G. Knepley   PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE));
3023157dc65SMatthew G. Knepley   PetscCall(VecSetUp(u));
3033157dc65SMatthew G. Knepley   PetscCall(TSSetSolution(ts, u));
3043157dc65SMatthew G. Knepley   PetscCall(VecDestroy(&u));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30640be0ff1SMatthew G. Knepley }
307db4653c2SDaniel Finn 
308d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetProblem(TS ts)
309d71ae5a4SJacob Faibussowitsch {
3103157dc65SMatthew G. Knepley   AppCtx *user;
3113157dc65SMatthew G. Knepley   DM      sw;
31240be0ff1SMatthew G. Knepley 
3133157dc65SMatthew G. Knepley   PetscFunctionBegin;
3143157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
3153157dc65SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, (void **)&user));
3163157dc65SMatthew G. Knepley   // Define unified system for (X, V)
3173157dc65SMatthew G. Knepley   {
3183157dc65SMatthew G. Knepley     Mat      J;
3193157dc65SMatthew G. Knepley     PetscInt dim, Np;
3203157dc65SMatthew G. Knepley 
3213157dc65SMatthew G. Knepley     PetscCall(DMGetDimension(sw, &dim));
3223157dc65SMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &Np));
3233157dc65SMatthew G. Knepley     PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
3243157dc65SMatthew G. Knepley     PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE));
3253157dc65SMatthew G. Knepley     PetscCall(MatSetBlockSize(J, 2 * dim));
3263157dc65SMatthew G. Knepley     PetscCall(MatSetFromOptions(J));
3273157dc65SMatthew G. Knepley     PetscCall(MatSetUp(J));
3283157dc65SMatthew G. Knepley     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user));
3293157dc65SMatthew G. Knepley     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user));
3303157dc65SMatthew G. Knepley     PetscCall(MatDestroy(&J));
33140be0ff1SMatthew G. Knepley   }
3323157dc65SMatthew G. Knepley   // Define split system for X and V
3333157dc65SMatthew G. Knepley   {
3343157dc65SMatthew G. Knepley     IS              isx, isv, istmp;
3353157dc65SMatthew G. Knepley     const PetscInt *idx;
3363157dc65SMatthew G. Knepley     PetscInt        dim, Np;
3373157dc65SMatthew G. Knepley 
3383157dc65SMatthew G. Knepley     PetscCall(DMGetDimension(sw, &dim));
3393157dc65SMatthew G. Knepley     PetscCall(DMSwarmGetLocalSize(sw, &Np));
3403157dc65SMatthew G. Knepley     PetscCall(ISCreateStride(PETSC_COMM_SELF, Np, 0, 2, &istmp));
3413157dc65SMatthew G. Knepley     PetscCall(ISGetIndices(istmp, &idx));
3423157dc65SMatthew G. Knepley     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, Np, idx, PETSC_COPY_VALUES, &isx));
3433157dc65SMatthew G. Knepley     PetscCall(ISRestoreIndices(istmp, &idx));
3443157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&istmp));
3453157dc65SMatthew G. Knepley     PetscCall(ISCreateStride(PETSC_COMM_SELF, Np, 1, 2, &istmp));
3463157dc65SMatthew G. Knepley     PetscCall(ISGetIndices(istmp, &idx));
3473157dc65SMatthew G. Knepley     PetscCall(ISCreateBlock(PETSC_COMM_SELF, dim, Np, idx, PETSC_COPY_VALUES, &isv));
3483157dc65SMatthew G. Knepley     PetscCall(ISRestoreIndices(istmp, &idx));
3493157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&istmp));
3503157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetIS(ts, "position", isx));
3513157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetIS(ts, "momentum", isv));
3523157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&isx));
3533157dc65SMatthew G. Knepley     PetscCall(ISDestroy(&isv));
3543157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user));
3553157dc65SMatthew G. Knepley     PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user));
3563157dc65SMatthew G. Knepley   }
3573157dc65SMatthew G. Knepley   // Define symplectic formulation U_t = S . G, where G = grad F
358d71ae5a4SJacob Faibussowitsch   {
359d71ae5a4SJacob Faibussowitsch     PetscCall(TSDiscGradSetFormulation(ts, RHSJacobianS, RHSObjectiveF, RHSFunctionG, user));
360d71ae5a4SJacob Faibussowitsch   }
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36240be0ff1SMatthew G. Knepley }
36340be0ff1SMatthew G. Knepley 
364d71ae5a4SJacob Faibussowitsch static PetscErrorCode InitializeSolve(TS ts, Vec u)
365d71ae5a4SJacob Faibussowitsch {
3663157dc65SMatthew G. Knepley   DM      sw;
3673157dc65SMatthew G. Knepley   Vec     gc, gc0;
3683157dc65SMatthew G. Knepley   IS      isx, isv;
3693157dc65SMatthew G. Knepley   AppCtx *user;
37040be0ff1SMatthew G. Knepley 
37140be0ff1SMatthew G. Knepley   PetscFunctionBeginUser;
3723157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
3733157dc65SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, &user));
3743157dc65SMatthew G. Knepley   {
3753157dc65SMatthew G. Knepley     PetscReal v0[1] = {1.};
376db4653c2SDaniel Finn 
3773157dc65SMatthew G. Knepley     PetscCall(DMSwarmInitializeCoordinates(sw));
3783157dc65SMatthew G. Knepley     PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0));
3793157dc65SMatthew G. Knepley     PetscCall(SetProblem(ts));
3803157dc65SMatthew G. Knepley   }
3813157dc65SMatthew G. Knepley   PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
3823157dc65SMatthew G. Knepley   PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
3833157dc65SMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
3843157dc65SMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0));
3853157dc65SMatthew G. Knepley   PetscCall(VecCopy(gc, gc0));
3863157dc65SMatthew G. Knepley   PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc));
3873157dc65SMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
3883157dc65SMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0));
3893157dc65SMatthew G. Knepley   PetscCall(VecISSet(u, isv, 0.));
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3913157dc65SMatthew G. Knepley }
3923157dc65SMatthew G. Knepley 
393d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
394d71ae5a4SJacob Faibussowitsch {
3953157dc65SMatthew G. Knepley   MPI_Comm           comm;
3963157dc65SMatthew G. Knepley   DM                 sw;
3973157dc65SMatthew G. Knepley   AppCtx            *user;
3983157dc65SMatthew G. Knepley   const PetscScalar *u;
3993157dc65SMatthew G. Knepley   const PetscReal   *coords;
4003157dc65SMatthew G. Knepley   PetscScalar       *e;
4013157dc65SMatthew G. Knepley   PetscReal          t;
4023157dc65SMatthew G. Knepley   PetscInt           dim, d, Np, p;
4033157dc65SMatthew G. Knepley 
4043157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
4053157dc65SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
4063157dc65SMatthew G. Knepley   PetscCall(TSGetDM(ts, &sw));
4073157dc65SMatthew G. Knepley   PetscCall(DMGetApplicationContext(sw, &user));
4083157dc65SMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
4093157dc65SMatthew G. Knepley   PetscCall(TSGetSolveTime(ts, &t));
4103157dc65SMatthew G. Knepley   PetscCall(VecGetArray(E, &e));
4119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
4129566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(U, &Np));
4133157dc65SMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
4143157dc65SMatthew G. Knepley   Np /= 2 * dim;
41540be0ff1SMatthew G. Knepley   for (p = 0; p < Np; ++p) {
4163157dc65SMatthew G. Knepley     const PetscReal omega = user->omega;
4173157dc65SMatthew G. Knepley     const PetscReal ct    = PetscCosReal(omega * t);
4183157dc65SMatthew G. Knepley     const PetscReal st    = PetscSinReal(omega * t);
4193157dc65SMatthew G. Knepley     const PetscReal x0    = DMPlex_NormD_Internal(dim, &coords[p * dim]);
4203157dc65SMatthew G. Knepley     const PetscReal ex    = x0 * ct;
4213157dc65SMatthew G. Knepley     const PetscReal ev    = -x0 * omega * st;
4223157dc65SMatthew G. Knepley     const PetscReal x     = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 0) * dim], &coords[p * dim]) / x0;
4233157dc65SMatthew G. Knepley     const PetscReal v     = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &coords[p * dim]) / x0;
4243157dc65SMatthew G. Knepley 
4253157dc65SMatthew G. Knepley     if (user->error) {
4263157dc65SMatthew G. Knepley       const PetscReal en   = 0.5 * (v * v + PetscSqr(omega) * x * x);
4273157dc65SMatthew G. Knepley       const PetscReal exen = 0.5 * PetscSqr(omega * x0);
4283771a240SStefano Zampini       PetscCall(PetscPrintf(comm, "p%" PetscInt_FMT " error [%.2f %.2f] 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)));
4293157dc65SMatthew G. Knepley     }
4303157dc65SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
4313157dc65SMatthew G. Knepley       e[(p * 2 + 0) * dim + d] = u[(p * 2 + 0) * dim + d] - coords[p * dim + d] * ct;
4323157dc65SMatthew G. Knepley       e[(p * 2 + 1) * dim + d] = u[(p * 2 + 1) * dim + d] + coords[p * dim + d] * omega * st;
4333157dc65SMatthew G. Knepley     }
4343157dc65SMatthew G. Knepley   }
4353157dc65SMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
4363157dc65SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(U, &u));
4373157dc65SMatthew G. Knepley   PetscCall(VecRestoreArray(E, &e));
4383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4393157dc65SMatthew G. Knepley }
4403157dc65SMatthew G. Knepley 
441d71ae5a4SJacob Faibussowitsch static PetscErrorCode EnergyMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
442d71ae5a4SJacob Faibussowitsch {
4433157dc65SMatthew G. Knepley   const PetscReal    omega = ((AppCtx *)ctx)->omega;
4443157dc65SMatthew G. Knepley   const PetscInt     ostep = ((AppCtx *)ctx)->ostep;
4453157dc65SMatthew G. Knepley   DM                 sw;
4463157dc65SMatthew G. Knepley   const PetscScalar *u;
4473157dc65SMatthew G. Knepley   PetscReal          dt;
4483157dc65SMatthew G. Knepley   PetscInt           dim, Np, p;
4493157dc65SMatthew G. Knepley   MPI_Comm           comm;
4503157dc65SMatthew G. Knepley 
4513157dc65SMatthew G. Knepley   PetscFunctionBeginUser;
4523157dc65SMatthew G. Knepley   if (step % ostep == 0) {
4533157dc65SMatthew G. Knepley     PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
4543157dc65SMatthew G. Knepley     PetscCall(TSGetDM(ts, &sw));
4553157dc65SMatthew G. Knepley     PetscCall(TSGetTimeStep(ts, &dt));
4563157dc65SMatthew G. Knepley     PetscCall(DMGetDimension(sw, &dim));
4573157dc65SMatthew G. Knepley     PetscCall(VecGetArrayRead(U, &u));
4583157dc65SMatthew G. Knepley     PetscCall(VecGetLocalSize(U, &Np));
4593157dc65SMatthew G. Knepley     Np /= 2 * dim;
4603157dc65SMatthew G. Knepley     if (!step) PetscCall(PetscPrintf(comm, "Time     Step Part     Energy Mod Energy\n"));
4613157dc65SMatthew G. Knepley     for (p = 0; p < Np; ++p) {
4623157dc65SMatthew G. Knepley       const PetscReal x2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 0) * dim], &u[(p * 2 + 0) * dim]);
4633157dc65SMatthew G. Knepley       const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]);
4643157dc65SMatthew G. Knepley       const PetscReal E  = 0.5 * (v2 + PetscSqr(omega) * x2);
4653157dc65SMatthew G. Knepley       const PetscReal mE = 0.5 * (v2 + PetscSqr(omega) * x2 - PetscSqr(omega) * dt * PetscSqrtReal(x2 * v2));
4663157dc65SMatthew G. Knepley 
46763a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "%.6lf %4" PetscInt_FMT " %4" PetscInt_FMT " %10.4lf %10.4lf\n", (double)t, step, p, (double)E, (double)mE));
46840be0ff1SMatthew G. Knepley     }
4699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U, &u));
4703157dc65SMatthew G. Knepley   }
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47240be0ff1SMatthew G. Knepley }
47340be0ff1SMatthew G. Knepley 
474d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
475d71ae5a4SJacob Faibussowitsch {
4763157dc65SMatthew G. Knepley   DM     dm, sw;
4773157dc65SMatthew G. Knepley   TS     ts;
4783157dc65SMatthew G. Knepley   Vec    u;
479c4762a1bSJed Brown   AppCtx user;
480c4762a1bSJed Brown 
481327415f7SBarry Smith   PetscFunctionBeginUser;
4829566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4833157dc65SMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
4843157dc65SMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
4853157dc65SMatthew G. Knepley   PetscCall(CreateSwarm(dm, &user, &sw));
4869566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(sw, &user));
487c4762a1bSJed Brown 
4883157dc65SMatthew G. Knepley   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
4899566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
4909566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, sw));
4919566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 0.1));
4929566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.00001));
4939566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 100));
4949566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
4953157dc65SMatthew G. Knepley   PetscCall(TSMonitorSet(ts, EnergyMonitor, &user, NULL));
4969566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
497*f940b0e3Sdanofinn 
4989566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
4999566063dSJacob Faibussowitsch   PetscCall(TSSetComputeExactError(ts, ComputeError));
50040be0ff1SMatthew G. Knepley 
5013157dc65SMatthew G. Knepley   PetscCall(CreateSolution(ts));
5023157dc65SMatthew G. Knepley   PetscCall(TSGetSolution(ts, &u));
5039566063dSJacob Faibussowitsch   PetscCall(TSComputeInitialCondition(ts, u));
5043157dc65SMatthew G. Knepley   PetscCall(TSSolve(ts, NULL));
505c4762a1bSJed Brown 
5069566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
5079566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
5089566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
5099566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
510b122ec5aSJacob Faibussowitsch   return 0;
511c4762a1bSJed Brown }
512c4762a1bSJed Brown 
513c4762a1bSJed Brown /*TEST
514c4762a1bSJed Brown 
515c4762a1bSJed Brown    build:
516c4762a1bSJed Brown      requires: triangle !single !complex
51740be0ff1SMatthew G. Knepley 
5183157dc65SMatthew G. Knepley    testset:
5193157dc65SMatthew G. Knepley      args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1 \
5203157dc65SMatthew G. Knepley            -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
5213157dc65SMatthew G. Knepley            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
5223157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
523c4762a1bSJed Brown      test:
5243157dc65SMatthew G. Knepley        suffix: bsi_1
5253157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 1
5263157dc65SMatthew G. Knepley      test:
5273157dc65SMatthew G. Knepley        suffix: bsi_2
5283157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 2
5293157dc65SMatthew G. Knepley      test:
5303157dc65SMatthew G. Knepley        suffix: bsi_3
5313157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 3
5323157dc65SMatthew G. Knepley      test:
5333157dc65SMatthew G. Knepley        suffix: bsi_4
5343157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
53540be0ff1SMatthew G. Knepley 
5363157dc65SMatthew G. Knepley    testset:
5373157dc65SMatthew 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 \
5383157dc65SMatthew G. Knepley            -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
5393157dc65SMatthew G. Knepley            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
5403157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
541c4762a1bSJed Brown      test:
5423157dc65SMatthew G. Knepley        suffix: bsi_2d_1
5433157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 1
5443157dc65SMatthew G. Knepley      test:
5453157dc65SMatthew G. Knepley        suffix: bsi_2d_2
5463157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 2
5473157dc65SMatthew G. Knepley      test:
5483157dc65SMatthew G. Knepley        suffix: bsi_2d_3
5493157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 3
5503157dc65SMatthew G. Knepley      test:
5513157dc65SMatthew G. Knepley        suffix: bsi_2d_4
5523157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
55340be0ff1SMatthew G. Knepley 
5543157dc65SMatthew G. Knepley    testset:
5553157dc65SMatthew 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 \
5563157dc65SMatthew G. Knepley            -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
5573157dc65SMatthew G. Knepley            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
5583157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
559c4762a1bSJed Brown      test:
5603157dc65SMatthew G. Knepley        suffix: bsi_3d_1
5613157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 1
5623157dc65SMatthew G. Knepley      test:
5633157dc65SMatthew G. Knepley        suffix: bsi_3d_2
5643157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 2
5653157dc65SMatthew G. Knepley      test:
5663157dc65SMatthew G. Knepley        suffix: bsi_3d_3
5673157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 3
5683157dc65SMatthew G. Knepley      test:
5693157dc65SMatthew G. Knepley        suffix: bsi_3d_4
5703157dc65SMatthew G. Knepley        args: -ts_basicsymplectic_type 4 -ts_dt 0.0001
57140be0ff1SMatthew G. Knepley 
5723157dc65SMatthew G. Knepley    testset:
5733157dc65SMatthew G. Knepley      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
5743157dc65SMatthew G. Knepley            -ts_type theta -ts_theta_theta 0.5 -ts_convergence_estimate -convest_num_refine 2 \
5753157dc65SMatthew G. Knepley              -mat_type baij -ksp_error_if_not_converged -pc_type lu \
5763157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
57740be0ff1SMatthew G. Knepley      test:
5783157dc65SMatthew G. Knepley        suffix: im_1d
5793157dc65SMatthew G. Knepley        args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1
5803157dc65SMatthew G. Knepley      test:
5813157dc65SMatthew G. Knepley        suffix: im_2d
5823157dc65SMatthew 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
5833157dc65SMatthew G. Knepley      test:
5843157dc65SMatthew G. Knepley        suffix: im_3d
5853157dc65SMatthew 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
58640be0ff1SMatthew G. Knepley 
5873157dc65SMatthew G. Knepley    testset:
5883157dc65SMatthew G. Knepley      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
589*f940b0e3Sdanofinn            -ts_type discgrad -ts_discgrad_type none -ts_convergence_estimate -convest_num_refine 2 \
5903157dc65SMatthew G. Knepley              -mat_type baij -ksp_error_if_not_converged -pc_type lu \
5913157dc65SMatthew G. Knepley            -dm_view -output_step 50 -error
59240be0ff1SMatthew G. Knepley      test:
5933157dc65SMatthew G. Knepley        suffix: dg_1d
5943157dc65SMatthew G. Knepley        args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1
595db4653c2SDaniel Finn      test:
5963157dc65SMatthew G. Knepley        suffix: dg_2d
5973157dc65SMatthew 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
5983157dc65SMatthew G. Knepley      test:
5993157dc65SMatthew G. Knepley        suffix: dg_3d
6003157dc65SMatthew 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
601c4762a1bSJed Brown 
602*f940b0e3Sdanofinn    testset:
603*f940b0e3Sdanofinn      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_density constant \
604*f940b0e3Sdanofinn             -ts_type discgrad -ts_convergence_estimate -convest_num_refine 2 \
605*f940b0e3Sdanofinn              -mat_type baij -ksp_error_if_not_converged -pc_type lu \
606*f940b0e3Sdanofinn             -dm_view -output_step 50 -error
607*f940b0e3Sdanofinn      test:
608*f940b0e3Sdanofinn        suffix: dg_gonzalez
609*f940b0e3Sdanofinn        args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1 -ts_discgrad_type gonzalez
610*f940b0e3Sdanofinn      test:
611*f940b0e3Sdanofinn        suffix: dg_average
612*f940b0e3Sdanofinn        args: -dm_plex_dim 1 -dm_plex_box_faces 1 -dm_plex_box_lower -1 -dm_plex_box_upper 1 -ts_discgrad_type average
613*f940b0e3Sdanofinn 
614c4762a1bSJed Brown TEST*/
615