xref: /petsc/src/ts/tests/ex80.c (revision 414ab8ebc8c11e9c1bafa4a96d65c2b23738fa87)
1*414ab8ebSDavid Kamensky static char help[] = "Constant acceleration check with 2nd-order generalized-alpha.\n";
2*414ab8ebSDavid Kamensky 
3*414ab8ebSDavid Kamensky #include <petscts.h>
4*414ab8ebSDavid Kamensky 
5*414ab8ebSDavid Kamensky typedef struct {
6*414ab8ebSDavid Kamensky   PetscReal a0;     /* contant acceleration  */
7*414ab8ebSDavid Kamensky   PetscReal u0, v0; /* initial conditions */
8*414ab8ebSDavid Kamensky   PetscReal radius; /* spectral radius of integrator */
9*414ab8ebSDavid Kamensky } UserParams;
10*414ab8ebSDavid Kamensky 
11*414ab8ebSDavid Kamensky static void Exact(PetscReal t, PetscReal a0, PetscReal u0, PetscReal v0, PetscScalar *ut, PetscScalar *vt)
12*414ab8ebSDavid Kamensky {
13*414ab8ebSDavid Kamensky   if (ut) *ut = u0 + v0 * t + 0.5 * a0 * t * t;
14*414ab8ebSDavid Kamensky   if (vt) *vt = v0 + a0 * t;
15*414ab8ebSDavid Kamensky }
16*414ab8ebSDavid Kamensky 
17*414ab8ebSDavid Kamensky PetscErrorCode Residual(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec R, void *ctx)
18*414ab8ebSDavid Kamensky {
19*414ab8ebSDavid Kamensky   UserParams        *user = (UserParams *)ctx;
20*414ab8ebSDavid Kamensky   const PetscScalar *a;
21*414ab8ebSDavid Kamensky   PetscScalar       *r;
22*414ab8ebSDavid Kamensky 
23*414ab8ebSDavid Kamensky   PetscFunctionBeginUser;
24*414ab8ebSDavid Kamensky   PetscCall(VecGetArrayRead(A, &a));
25*414ab8ebSDavid Kamensky   PetscCall(VecGetArrayWrite(R, &r));
26*414ab8ebSDavid Kamensky 
27*414ab8ebSDavid Kamensky   r[0] = a[0] - user->a0;
28*414ab8ebSDavid Kamensky 
29*414ab8ebSDavid Kamensky   PetscCall(VecRestoreArrayRead(A, &a));
30*414ab8ebSDavid Kamensky   PetscCall(VecRestoreArrayWrite(R, &r));
31*414ab8ebSDavid Kamensky   PetscFunctionReturn(PETSC_SUCCESS);
32*414ab8ebSDavid Kamensky }
33*414ab8ebSDavid Kamensky 
34*414ab8ebSDavid Kamensky PetscErrorCode Tangent(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P, void *ctx)
35*414ab8ebSDavid Kamensky {
36*414ab8ebSDavid Kamensky   PetscReal T = 0;
37*414ab8ebSDavid Kamensky 
38*414ab8ebSDavid Kamensky   PetscFunctionBeginUser;
39*414ab8ebSDavid Kamensky   T = shiftA;
40*414ab8ebSDavid Kamensky 
41*414ab8ebSDavid Kamensky   PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES));
42*414ab8ebSDavid Kamensky   PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
43*414ab8ebSDavid Kamensky   PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
44*414ab8ebSDavid Kamensky   if (J != P) {
45*414ab8ebSDavid Kamensky     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
46*414ab8ebSDavid Kamensky     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
47*414ab8ebSDavid Kamensky   }
48*414ab8ebSDavid Kamensky   PetscFunctionReturn(PETSC_SUCCESS);
49*414ab8ebSDavid Kamensky }
50*414ab8ebSDavid Kamensky 
51*414ab8ebSDavid Kamensky int main(int argc, char *argv[])
52*414ab8ebSDavid Kamensky {
53*414ab8ebSDavid Kamensky   PetscMPIInt  size;
54*414ab8ebSDavid Kamensky   TS           ts;
55*414ab8ebSDavid Kamensky   Vec          R;
56*414ab8ebSDavid Kamensky   Mat          J;
57*414ab8ebSDavid Kamensky   Vec          U, V;
58*414ab8ebSDavid Kamensky   PetscScalar *u, *v, u_exact, v_exact;
59*414ab8ebSDavid Kamensky   PetscReal    u_err, v_err;
60*414ab8ebSDavid Kamensky   PetscReal    atol    = 1e-15;
61*414ab8ebSDavid Kamensky   PetscReal    t_final = 3.0;
62*414ab8ebSDavid Kamensky   PetscInt     n_step  = 8;
63*414ab8ebSDavid Kamensky   UserParams   user    = {/*a0=*/1, /*u0=*/1, /*v0=*/0, /*radius=*/0.0};
64*414ab8ebSDavid Kamensky 
65*414ab8ebSDavid Kamensky   PetscFunctionBeginUser;
66*414ab8ebSDavid Kamensky   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
67*414ab8ebSDavid Kamensky   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
68*414ab8ebSDavid Kamensky   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
69*414ab8ebSDavid Kamensky 
70*414ab8ebSDavid Kamensky   PetscOptionsBegin(PETSC_COMM_SELF, "", "ex80 options", "");
71*414ab8ebSDavid Kamensky   PetscCall(PetscOptionsReal("-acceleration", "Constant acceleration", __FILE__, user.a0, &user.a0, NULL));
72*414ab8ebSDavid Kamensky   PetscCall(PetscOptionsReal("-initial_u", "Initial displacement", __FILE__, user.u0, &user.u0, NULL));
73*414ab8ebSDavid Kamensky   PetscCall(PetscOptionsReal("-initial_v", "Initial velocity", __FILE__, user.v0, &user.v0, NULL));
74*414ab8ebSDavid Kamensky   PetscCall(PetscOptionsReal("-radius", "Spectral radius", __FILE__, user.radius, &user.radius, NULL));
75*414ab8ebSDavid Kamensky   PetscOptionsEnd();
76*414ab8ebSDavid Kamensky 
77*414ab8ebSDavid Kamensky   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
78*414ab8ebSDavid Kamensky   PetscCall(TSSetType(ts, TSALPHA2));
79*414ab8ebSDavid Kamensky   PetscCall(TSSetMaxTime(ts, t_final));
80*414ab8ebSDavid Kamensky   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
81*414ab8ebSDavid Kamensky   PetscCall(TSSetTimeStep(ts, t_final / n_step));
82*414ab8ebSDavid Kamensky   PetscCall(TSAlpha2SetRadius(ts, user.radius));
83*414ab8ebSDavid Kamensky 
84*414ab8ebSDavid Kamensky   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &R));
85*414ab8ebSDavid Kamensky   PetscCall(VecSetUp(R));
86*414ab8ebSDavid Kamensky   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 1, 1, NULL, &J));
87*414ab8ebSDavid Kamensky   PetscCall(MatSetUp(J));
88*414ab8ebSDavid Kamensky   PetscCall(TSSetI2Function(ts, R, Residual, &user));
89*414ab8ebSDavid Kamensky   PetscCall(TSSetI2Jacobian(ts, J, J, Tangent, &user));
90*414ab8ebSDavid Kamensky   PetscCall(VecDestroy(&R));
91*414ab8ebSDavid Kamensky   PetscCall(MatDestroy(&J));
92*414ab8ebSDavid Kamensky 
93*414ab8ebSDavid Kamensky   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &U));
94*414ab8ebSDavid Kamensky   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &V));
95*414ab8ebSDavid Kamensky   PetscCall(VecGetArrayWrite(U, &u));
96*414ab8ebSDavid Kamensky   PetscCall(VecGetArrayWrite(V, &v));
97*414ab8ebSDavid Kamensky   u[0] = user.u0;
98*414ab8ebSDavid Kamensky   v[0] = user.v0;
99*414ab8ebSDavid Kamensky   PetscCall(VecRestoreArrayWrite(U, &u));
100*414ab8ebSDavid Kamensky   PetscCall(VecRestoreArrayWrite(V, &v));
101*414ab8ebSDavid Kamensky 
102*414ab8ebSDavid Kamensky   PetscCall(TS2SetSolution(ts, U, V));
103*414ab8ebSDavid Kamensky   PetscCall(TSSetFromOptions(ts));
104*414ab8ebSDavid Kamensky   PetscCall(TSSolve(ts, NULL));
105*414ab8ebSDavid Kamensky 
106*414ab8ebSDavid Kamensky   PetscCall(VecGetArray(U, &u));
107*414ab8ebSDavid Kamensky   PetscCall(VecGetArray(V, &v));
108*414ab8ebSDavid Kamensky   Exact(t_final, user.a0, user.u0, user.v0, &u_exact, &v_exact);
109*414ab8ebSDavid Kamensky   u_err = PetscAbsScalar(u[0] - u_exact);
110*414ab8ebSDavid Kamensky   v_err = PetscAbsScalar(v[0] - v_exact);
111*414ab8ebSDavid Kamensky   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "u(t=%g) = %g (error = %g)\n", (double)t_final, (double)PetscRealPart(u[0]), (double)u_err));
112*414ab8ebSDavid Kamensky   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "v(t=%g) = %g (error = %g)\n", (double)t_final, (double)PetscRealPart(v[0]), (double)v_err));
113*414ab8ebSDavid Kamensky   PetscCheck(u_err < atol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Inexact displacement.");
114*414ab8ebSDavid Kamensky   PetscCheck(v_err < atol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Inexact velocity.");
115*414ab8ebSDavid Kamensky   PetscCall(VecRestoreArray(U, &u));
116*414ab8ebSDavid Kamensky   PetscCall(VecRestoreArray(V, &v));
117*414ab8ebSDavid Kamensky 
118*414ab8ebSDavid Kamensky   PetscCall(VecDestroy(&U));
119*414ab8ebSDavid Kamensky   PetscCall(VecDestroy(&V));
120*414ab8ebSDavid Kamensky   PetscCall(TSDestroy(&ts));
121*414ab8ebSDavid Kamensky   PetscCall(PetscFinalize());
122*414ab8ebSDavid Kamensky   return 0;
123*414ab8ebSDavid Kamensky }
124*414ab8ebSDavid Kamensky 
125*414ab8ebSDavid Kamensky /*TEST
126*414ab8ebSDavid Kamensky 
127*414ab8ebSDavid Kamensky     test:
128*414ab8ebSDavid Kamensky       suffix: a
129*414ab8ebSDavid Kamensky       args: -radius 0.0
130*414ab8ebSDavid Kamensky       requires: !single
131*414ab8ebSDavid Kamensky 
132*414ab8ebSDavid Kamensky     test:
133*414ab8ebSDavid Kamensky       suffix: b
134*414ab8ebSDavid Kamensky       args: -radius 0.5
135*414ab8ebSDavid Kamensky       requires: !single
136*414ab8ebSDavid Kamensky 
137*414ab8ebSDavid Kamensky     test:
138*414ab8ebSDavid Kamensky       suffix: c
139*414ab8ebSDavid Kamensky       args: -radius 1.0
140*414ab8ebSDavid Kamensky       requires: !single
141*414ab8ebSDavid Kamensky 
142*414ab8ebSDavid Kamensky TEST*/
143