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