xref: /petsc/src/ts/tutorials/ex43.c (revision 220f924abfb7c127e3e73e9ee375e6c993b7366b)
1c4762a1bSJed Brown static char help[] = "Single-DOF oscillator formulated as a second-order system.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscts.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   PetscReal Omega;    /* natural frequency */
7c4762a1bSJed Brown   PetscReal Xi;       /* damping coefficient  */
8c4762a1bSJed Brown   PetscReal u0, v0;   /* initial conditions */
9*220f924aSDavid Kamensky   PetscBool use_pred; /* whether to use a predictor callback */
10c4762a1bSJed Brown } UserParams;
11c4762a1bSJed Brown 
12d71ae5a4SJacob Faibussowitsch static void Exact(PetscReal t, PetscReal omega, PetscReal xi, PetscReal u0, PetscReal v0, PetscReal *ut, PetscReal *vt)
13d71ae5a4SJacob Faibussowitsch {
14c4762a1bSJed Brown   PetscReal u, v;
15c4762a1bSJed Brown   if (xi < 1) {
16c4762a1bSJed Brown     PetscReal a  = xi * omega;
17303a5415SBarry Smith     PetscReal w  = PetscSqrtReal(1 - xi * xi) * omega;
18c4762a1bSJed Brown     PetscReal C1 = (v0 + a * u0) / w;
19c4762a1bSJed Brown     PetscReal C2 = u0;
20303a5415SBarry Smith     u            = PetscExpReal(-a * t) * (C1 * PetscSinReal(w * t) + C2 * PetscCosReal(w * t));
21303a5415SBarry Smith     v            = (-a * PetscExpReal(-a * t) * (C1 * PetscSinReal(w * t) + C2 * PetscCosReal(w * t)) + w * PetscExpReal(-a * t) * (C1 * PetscCosReal(w * t) - C2 * PetscSinReal(w * t)));
22c4762a1bSJed Brown   } else if (xi > 1) {
23303a5415SBarry Smith     PetscReal w  = PetscSqrtReal(xi * xi - 1) * omega;
24c4762a1bSJed Brown     PetscReal C1 = (w * u0 + xi * u0 + v0) / (2 * w);
25c4762a1bSJed Brown     PetscReal C2 = (w * u0 - xi * u0 - v0) / (2 * w);
26303a5415SBarry Smith     u            = C1 * PetscExpReal((-xi + w) * t) + C2 * PetscExpReal((-xi - w) * t);
27303a5415SBarry Smith     v            = C1 * (-xi + w) * PetscExpReal((-xi + w) * t) + C2 * (-xi - w) * PetscExpReal((-xi - w) * t);
28c4762a1bSJed Brown   } else {
29c4762a1bSJed Brown     PetscReal a  = xi * omega;
30c4762a1bSJed Brown     PetscReal C1 = v0 + a * u0;
31c4762a1bSJed Brown     PetscReal C2 = u0;
32303a5415SBarry Smith     u            = (C1 * t + C2) * PetscExpReal(-a * t);
33303a5415SBarry Smith     v            = (C1 - a * (C1 * t + C2)) * PetscExpReal(-a * t);
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown   if (ut) *ut = u;
36c4762a1bSJed Brown   if (vt) *vt = v;
37c4762a1bSJed Brown }
38c4762a1bSJed Brown 
39d71ae5a4SJacob Faibussowitsch PetscErrorCode Solution(TS ts, PetscReal t, Vec X, void *ctx)
40d71ae5a4SJacob Faibussowitsch {
41c4762a1bSJed Brown   UserParams  *user = (UserParams *)ctx;
42c4762a1bSJed Brown   PetscReal    u, v;
43c4762a1bSJed Brown   PetscScalar *x;
44c4762a1bSJed Brown 
457510d9b0SBarry Smith   PetscFunctionBeginUser;
46c4762a1bSJed Brown   Exact(t, user->Omega, user->Xi, user->u0, user->v0, &u, &v);
479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
48c4762a1bSJed Brown   x[0] = (PetscScalar)u;
499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53*220f924aSDavid Kamensky PetscErrorCode Predictor(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx)
54*220f924aSDavid Kamensky {
55*220f924aSDavid Kamensky   PetscReal dt, accel_fac;
56*220f924aSDavid Kamensky 
57*220f924aSDavid Kamensky   PetscFunctionBeginUser;
58*220f924aSDavid Kamensky   PetscCall(TSGetTimeStep(ts, &dt));
59*220f924aSDavid Kamensky   accel_fac = 0.5 * dt * dt;
60*220f924aSDavid Kamensky   /* X1 = X0 + dt*V0 + accel_fac*A0 */
61*220f924aSDavid Kamensky   PetscCall(VecCopy(X0, X1));
62*220f924aSDavid Kamensky   PetscCall(VecAXPY(X1, dt, V0));
63*220f924aSDavid Kamensky   PetscCall(VecAXPY(X1, accel_fac, A0));
64*220f924aSDavid Kamensky   PetscFunctionReturn(PETSC_SUCCESS);
65*220f924aSDavid Kamensky }
66*220f924aSDavid Kamensky 
67d71ae5a4SJacob Faibussowitsch PetscErrorCode Residual1(TS ts, PetscReal t, Vec U, Vec A, Vec R, void *ctx)
68d71ae5a4SJacob Faibussowitsch {
69c4762a1bSJed Brown   UserParams        *user  = (UserParams *)ctx;
70c4762a1bSJed Brown   PetscReal          Omega = user->Omega;
71c4762a1bSJed Brown   const PetscScalar *u, *a;
72c4762a1bSJed Brown   PetscScalar       *r;
73c4762a1bSJed Brown 
747510d9b0SBarry Smith   PetscFunctionBeginUser;
759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A, &a));
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(R, &r));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   r[0] = a[0] + (Omega * Omega) * u[0];
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A, &a));
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(R, &r));
849566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(R));
859566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(R));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89d71ae5a4SJacob Faibussowitsch PetscErrorCode Tangent1(TS ts, PetscReal t, Vec U, Vec A, PetscReal shiftA, Mat J, Mat P, void *ctx)
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown   UserParams *user  = (UserParams *)ctx;
92c4762a1bSJed Brown   PetscReal   Omega = user->Omega;
93c4762a1bSJed Brown   PetscReal   T     = 0;
94c4762a1bSJed Brown 
957510d9b0SBarry Smith   PetscFunctionBeginUser;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   T = shiftA + (Omega * Omega);
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES));
1009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
1019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
102c4762a1bSJed Brown   if (J != P) {
1039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1049566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
105c4762a1bSJed Brown   }
1063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109d71ae5a4SJacob Faibussowitsch PetscErrorCode Residual2(TS ts, PetscReal t, Vec U, Vec V, Vec A, Vec R, void *ctx)
110d71ae5a4SJacob Faibussowitsch {
111c4762a1bSJed Brown   UserParams        *user  = (UserParams *)ctx;
112c4762a1bSJed Brown   PetscReal          Omega = user->Omega, Xi = user->Xi;
113c4762a1bSJed Brown   const PetscScalar *u, *v, *a;
114c4762a1bSJed Brown   PetscScalar       *r;
115c4762a1bSJed Brown 
1167510d9b0SBarry Smith   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V, &v));
1199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A, &a));
1209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(R, &r));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   r[0] = a[0] + (2 * Xi * Omega) * v[0] + (Omega * Omega) * u[0];
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
1269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A, &a));
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(R, &r));
1289566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(R));
1299566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(R));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133d71ae5a4SJacob Faibussowitsch PetscErrorCode Tangent2(TS ts, PetscReal t, Vec U, Vec V, Vec A, PetscReal shiftV, PetscReal shiftA, Mat J, Mat P, void *ctx)
134d71ae5a4SJacob Faibussowitsch {
135c4762a1bSJed Brown   UserParams *user  = (UserParams *)ctx;
136c4762a1bSJed Brown   PetscReal   Omega = user->Omega, Xi = user->Xi;
137c4762a1bSJed Brown   PetscReal   T = 0;
138c4762a1bSJed Brown 
1397510d9b0SBarry Smith   PetscFunctionBeginUser;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   T = shiftA + shiftV * (2 * Xi * Omega) + (Omega * Omega);
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 0, 0, T, INSERT_VALUES));
1449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
1459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
146c4762a1bSJed Brown   if (J != P) {
1479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1489566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
149c4762a1bSJed Brown   }
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown   PetscMPIInt  size;
156c4762a1bSJed Brown   TS           ts;
157c4762a1bSJed Brown   Vec          R;
158c4762a1bSJed Brown   Mat          J;
159c4762a1bSJed Brown   Vec          U, V;
160c4762a1bSJed Brown   PetscScalar *u, *v;
161*220f924aSDavid Kamensky   UserParams   user = {/*Omega=*/1, /*Xi=*/0, /*u0=*/1, /*,v0=*/0, /*,use_pred=*/PETSC_FALSE};
162c4762a1bSJed Brown 
163327415f7SBarry Smith   PetscFunctionBeginUser;
1649566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1663c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
167c4762a1bSJed Brown 
168d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "ex43 options", "");
1699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-frequency", "Natual frequency", __FILE__, user.Omega, &user.Omega, NULL));
1709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-damping", "Damping coefficient", __FILE__, user.Xi, &user.Xi, NULL));
1719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-initial_u", "Initial displacement", __FILE__, user.u0, &user.u0, NULL));
1729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-initial_v", "Initial velocity", __FILE__, user.v0, &user.v0, NULL));
173*220f924aSDavid Kamensky   PetscCall(PetscOptionsBool("-use_pred", "Use a custom predictor", __FILE__, user.use_pred, &user.use_pred, NULL));
174d0609cedSBarry Smith   PetscOptionsEnd();
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1779566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSALPHA2));
1789566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 5 * (2 * PETSC_PI)));
1799566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1809566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.01));
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &R));
1839566063dSJacob Faibussowitsch   PetscCall(VecSetUp(R));
1849566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, 1, 1, NULL, &J));
1859566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
186c4762a1bSJed Brown   if (user.Xi) {
1879566063dSJacob Faibussowitsch     PetscCall(TSSetI2Function(ts, R, Residual2, &user));
1889566063dSJacob Faibussowitsch     PetscCall(TSSetI2Jacobian(ts, J, J, Tangent2, &user));
189c4762a1bSJed Brown   } else {
1909566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, R, Residual1, &user));
1919566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, J, J, Tangent1, &user));
192c4762a1bSJed Brown   }
1939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
1949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts, Solution, &user));
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &U));
1989566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &V));
1999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(U, &u));
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(V, &v));
201c4762a1bSJed Brown   u[0] = user.u0;
202c4762a1bSJed Brown   v[0] = user.v0;
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(U, &u));
2049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(V, &v));
205c4762a1bSJed Brown 
206*220f924aSDavid Kamensky   if (user.use_pred) PetscCall(TSAlpha2SetPredictor(ts, Predictor, NULL));
207*220f924aSDavid Kamensky 
2089566063dSJacob Faibussowitsch   PetscCall(TS2SetSolution(ts, U, V));
2099566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2109566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, NULL));
211c4762a1bSJed Brown 
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V));
2149566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2159566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
216b122ec5aSJacob Faibussowitsch   return 0;
217c4762a1bSJed Brown }
218c4762a1bSJed Brown 
219c4762a1bSJed Brown /*TEST
220c4762a1bSJed Brown 
221c4762a1bSJed Brown     test:
222c4762a1bSJed Brown       suffix: a
223*220f924aSDavid Kamensky       args: -ts_max_steps 10 -ts_view -use_pred
224c4762a1bSJed Brown       requires: !single
225c4762a1bSJed Brown 
226c4762a1bSJed Brown     test:
227c4762a1bSJed Brown       suffix: b
228c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_rtol 0 -ts_atol 1e-5 -ts_adapt_type basic -ts_adapt_monitor
229c4762a1bSJed Brown       requires: !single
230c4762a1bSJed Brown 
231c4762a1bSJed Brown TEST*/
232