1414ab8ebSDavid Kamensky static char help[] = "Constant acceleration check with 2nd-order generalized-alpha.\n";
2414ab8ebSDavid Kamensky
3414ab8ebSDavid Kamensky #include <petscts.h>
4414ab8ebSDavid Kamensky
5414ab8ebSDavid Kamensky typedef struct {
6d7c1f440SPierre Jolivet PetscReal a0; /* constant acceleration */
7414ab8ebSDavid Kamensky PetscReal u0, v0; /* initial conditions */
8414ab8ebSDavid Kamensky PetscReal radius; /* spectral radius of integrator */
9414ab8ebSDavid Kamensky } UserParams;
10414ab8ebSDavid Kamensky
Exact(PetscReal t,PetscReal a0,PetscReal u0,PetscReal v0,PetscScalar * ut,PetscScalar * vt)11414ab8ebSDavid Kamensky static void Exact(PetscReal t, PetscReal a0, PetscReal u0, PetscReal v0, PetscScalar *ut, PetscScalar *vt)
12414ab8ebSDavid Kamensky {
13414ab8ebSDavid Kamensky if (ut) *ut = u0 + v0 * t + 0.5 * a0 * t * t;
14414ab8ebSDavid Kamensky if (vt) *vt = v0 + a0 * t;
15414ab8ebSDavid Kamensky }
16414ab8ebSDavid Kamensky
Residual(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec R,PetscCtx ctx)17*2a8381b2SBarry Smith PetscErrorCode Residual(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec R, PetscCtx ctx)
18414ab8ebSDavid Kamensky {
19414ab8ebSDavid Kamensky UserParams *user = (UserParams *)ctx;
20414ab8ebSDavid Kamensky const PetscScalar *a;
21414ab8ebSDavid Kamensky PetscScalar *r;
22414ab8ebSDavid Kamensky
23414ab8ebSDavid Kamensky PetscFunctionBeginUser;
24414ab8ebSDavid Kamensky PetscCall(VecGetArrayRead(A, &a));
25414ab8ebSDavid Kamensky PetscCall(VecGetArrayWrite(R, &r));
26414ab8ebSDavid Kamensky
27414ab8ebSDavid Kamensky r[0] = a[0] - user->a0;
28414ab8ebSDavid Kamensky
29414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayRead(A, &a));
30414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayWrite(R, &r));
31414ab8ebSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS);
32414ab8ebSDavid Kamensky }
33414ab8ebSDavid Kamensky
Tangent(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P,PetscCtx ctx)34*2a8381b2SBarry Smith PetscErrorCode Tangent(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P, PetscCtx ctx)
35414ab8ebSDavid Kamensky {
36414ab8ebSDavid Kamensky PetscReal T = 0;
37414ab8ebSDavid Kamensky
38414ab8ebSDavid Kamensky PetscFunctionBeginUser;
39414ab8ebSDavid Kamensky T = shiftA;
40414ab8ebSDavid Kamensky
41414ab8ebSDavid Kamensky PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES));
42414ab8ebSDavid Kamensky PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
43414ab8ebSDavid Kamensky PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
44414ab8ebSDavid Kamensky if (J != P) {
45414ab8ebSDavid Kamensky PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
46414ab8ebSDavid Kamensky PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
47414ab8ebSDavid Kamensky }
48414ab8ebSDavid Kamensky PetscFunctionReturn(PETSC_SUCCESS);
49414ab8ebSDavid Kamensky }
50414ab8ebSDavid Kamensky
main(int argc,char * argv[])51414ab8ebSDavid Kamensky int main(int argc, char *argv[])
52414ab8ebSDavid Kamensky {
53414ab8ebSDavid Kamensky PetscMPIInt size;
54414ab8ebSDavid Kamensky TS ts;
55414ab8ebSDavid Kamensky Vec R;
56414ab8ebSDavid Kamensky Mat J;
57414ab8ebSDavid Kamensky Vec U, V;
58414ab8ebSDavid Kamensky PetscScalar *u, *v, u_exact, v_exact;
59414ab8ebSDavid Kamensky PetscReal u_err, v_err;
60414ab8ebSDavid Kamensky PetscReal atol = 1e-15;
61414ab8ebSDavid Kamensky PetscReal t_final = 3.0;
62414ab8ebSDavid Kamensky PetscInt n_step = 8;
63414ab8ebSDavid Kamensky UserParams user = {/*a0=*/1, /*u0=*/1, /*v0=*/0, /*radius=*/0.0};
64414ab8ebSDavid Kamensky
65414ab8ebSDavid Kamensky PetscFunctionBeginUser;
66414ab8ebSDavid Kamensky PetscCall(PetscInitialize(&argc, &argv, NULL, help));
67414ab8ebSDavid Kamensky PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
68414ab8ebSDavid Kamensky PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
69414ab8ebSDavid Kamensky
70414ab8ebSDavid Kamensky PetscOptionsBegin(PETSC_COMM_SELF, "", "ex80 options", "");
71414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-acceleration", "Constant acceleration", __FILE__, user.a0, &user.a0, NULL));
72414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-initial_u", "Initial displacement", __FILE__, user.u0, &user.u0, NULL));
73414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-initial_v", "Initial velocity", __FILE__, user.v0, &user.v0, NULL));
74414ab8ebSDavid Kamensky PetscCall(PetscOptionsReal("-radius", "Spectral radius", __FILE__, user.radius, &user.radius, NULL));
75414ab8ebSDavid Kamensky PetscOptionsEnd();
76414ab8ebSDavid Kamensky
77414ab8ebSDavid Kamensky PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
78414ab8ebSDavid Kamensky PetscCall(TSSetType(ts, TSALPHA2));
79414ab8ebSDavid Kamensky PetscCall(TSSetMaxTime(ts, t_final));
80414ab8ebSDavid Kamensky PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
81414ab8ebSDavid Kamensky PetscCall(TSSetTimeStep(ts, t_final / n_step));
82414ab8ebSDavid Kamensky PetscCall(TSAlpha2SetRadius(ts, user.radius));
83414ab8ebSDavid Kamensky
84414ab8ebSDavid Kamensky PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &R));
85414ab8ebSDavid Kamensky PetscCall(VecSetUp(R));
86414ab8ebSDavid Kamensky PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 1, 1, NULL, &J));
87414ab8ebSDavid Kamensky PetscCall(MatSetUp(J));
88414ab8ebSDavid Kamensky PetscCall(TSSetI2Function(ts, R, Residual, &user));
89414ab8ebSDavid Kamensky PetscCall(TSSetI2Jacobian(ts, J, J, Tangent, &user));
90414ab8ebSDavid Kamensky PetscCall(VecDestroy(&R));
91414ab8ebSDavid Kamensky PetscCall(MatDestroy(&J));
92414ab8ebSDavid Kamensky
93414ab8ebSDavid Kamensky PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &U));
94414ab8ebSDavid Kamensky PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &V));
95414ab8ebSDavid Kamensky PetscCall(VecGetArrayWrite(U, &u));
96414ab8ebSDavid Kamensky PetscCall(VecGetArrayWrite(V, &v));
97414ab8ebSDavid Kamensky u[0] = user.u0;
98414ab8ebSDavid Kamensky v[0] = user.v0;
99414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayWrite(U, &u));
100414ab8ebSDavid Kamensky PetscCall(VecRestoreArrayWrite(V, &v));
101414ab8ebSDavid Kamensky
102414ab8ebSDavid Kamensky PetscCall(TS2SetSolution(ts, U, V));
103414ab8ebSDavid Kamensky PetscCall(TSSetFromOptions(ts));
104414ab8ebSDavid Kamensky PetscCall(TSSolve(ts, NULL));
105414ab8ebSDavid Kamensky
106414ab8ebSDavid Kamensky PetscCall(VecGetArray(U, &u));
107414ab8ebSDavid Kamensky PetscCall(VecGetArray(V, &v));
108414ab8ebSDavid Kamensky Exact(t_final, user.a0, user.u0, user.v0, &u_exact, &v_exact);
109414ab8ebSDavid Kamensky u_err = PetscAbsScalar(u[0] - u_exact);
110414ab8ebSDavid Kamensky v_err = PetscAbsScalar(v[0] - v_exact);
111414ab8ebSDavid Kamensky PetscCall(PetscPrintf(PETSC_COMM_WORLD, "u(t=%g) = %g (error = %g)\n", (double)t_final, (double)PetscRealPart(u[0]), (double)u_err));
112414ab8ebSDavid Kamensky PetscCall(PetscPrintf(PETSC_COMM_WORLD, "v(t=%g) = %g (error = %g)\n", (double)t_final, (double)PetscRealPart(v[0]), (double)v_err));
113414ab8ebSDavid Kamensky PetscCheck(u_err < atol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Inexact displacement.");
114414ab8ebSDavid Kamensky PetscCheck(v_err < atol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Inexact velocity.");
115414ab8ebSDavid Kamensky PetscCall(VecRestoreArray(U, &u));
116414ab8ebSDavid Kamensky PetscCall(VecRestoreArray(V, &v));
117414ab8ebSDavid Kamensky
118414ab8ebSDavid Kamensky PetscCall(VecDestroy(&U));
119414ab8ebSDavid Kamensky PetscCall(VecDestroy(&V));
120414ab8ebSDavid Kamensky PetscCall(TSDestroy(&ts));
121414ab8ebSDavid Kamensky PetscCall(PetscFinalize());
122414ab8ebSDavid Kamensky return 0;
123414ab8ebSDavid Kamensky }
124414ab8ebSDavid Kamensky
125414ab8ebSDavid Kamensky /*TEST
126414ab8ebSDavid Kamensky
127414ab8ebSDavid Kamensky test:
128414ab8ebSDavid Kamensky suffix: a
129414ab8ebSDavid Kamensky args: -radius 0.0
130414ab8ebSDavid Kamensky requires: !single
131414ab8ebSDavid Kamensky
132414ab8ebSDavid Kamensky test:
133414ab8ebSDavid Kamensky suffix: b
134414ab8ebSDavid Kamensky args: -radius 0.5
135414ab8ebSDavid Kamensky requires: !single
136414ab8ebSDavid Kamensky
137414ab8ebSDavid Kamensky test:
138414ab8ebSDavid Kamensky suffix: c
139414ab8ebSDavid Kamensky args: -radius 1.0
140414ab8ebSDavid Kamensky requires: !single
141414ab8ebSDavid Kamensky
142414ab8ebSDavid Kamensky TEST*/
143