xref: /petsc/src/ts/tutorials/ex43.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 */
9c4762a1bSJed Brown } UserParams;
10c4762a1bSJed Brown 
11303a5415SBarry Smith static void Exact(PetscReal t,PetscReal omega,PetscReal xi,PetscReal u0,PetscReal v0,PetscReal *ut,PetscReal *vt)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscReal u,v;
14c4762a1bSJed Brown   if (xi < 1) {
15c4762a1bSJed Brown     PetscReal a  = xi*omega;
16303a5415SBarry Smith     PetscReal w  = PetscSqrtReal(1-xi*xi)*omega;
17c4762a1bSJed Brown     PetscReal C1 = (v0 + a*u0)/w;
18c4762a1bSJed Brown     PetscReal C2 = u0;
19303a5415SBarry Smith     u = PetscExpReal(-a*t) * (C1*PetscSinReal(w*t) + C2*PetscCosReal(w*t));
20303a5415SBarry 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)));
21c4762a1bSJed Brown   } else if (xi > 1) {
22303a5415SBarry Smith     PetscReal w  = PetscSqrtReal(xi*xi-1)*omega;
23c4762a1bSJed Brown     PetscReal C1 = (w*u0 + xi*u0 + v0)/(2*w);
24c4762a1bSJed Brown     PetscReal C2 = (w*u0 - xi*u0 - v0)/(2*w);
25303a5415SBarry Smith     u = C1*PetscExpReal((-xi+w)*t) + C2*PetscExpReal((-xi-w)*t);
26303a5415SBarry Smith     v = C1*(-xi+w)*PetscExpReal((-xi+w)*t) + C2*(-xi-w)*PetscExpReal((-xi-w)*t);
27c4762a1bSJed Brown   } else {
28c4762a1bSJed Brown     PetscReal a  = xi*omega;
29c4762a1bSJed Brown     PetscReal C1 = v0 + a*u0;
30c4762a1bSJed Brown     PetscReal C2 = u0;
31303a5415SBarry Smith     u = (C1*t + C2) * PetscExpReal(-a*t);
32303a5415SBarry Smith     v = (C1 - a*(C1*t + C2)) * PetscExpReal(-a*t);
33c4762a1bSJed Brown   }
34c4762a1bSJed Brown   if (ut) *ut = u;
35c4762a1bSJed Brown   if (vt) *vt = v;
36c4762a1bSJed Brown }
37c4762a1bSJed Brown 
38c4762a1bSJed Brown PetscErrorCode Solution(TS ts,PetscReal t,Vec X,void *ctx)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
41c4762a1bSJed Brown   PetscReal      u,v;
42c4762a1bSJed Brown   PetscScalar    *x;
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   PetscFunctionBegin;
45c4762a1bSJed Brown   Exact(t,user->Omega,user->Xi,user->u0,user->v0,&u,&v);
469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
47c4762a1bSJed Brown   x[0] = (PetscScalar)u;
489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
49c4762a1bSJed Brown   PetscFunctionReturn(0);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown PetscErrorCode Residual1(TS ts,PetscReal t,Vec U,Vec A,Vec R,void *ctx)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   UserParams        *user = (UserParams*)ctx;
55c4762a1bSJed Brown   PetscReal         Omega = user->Omega;
56c4762a1bSJed Brown   const PetscScalar *u,*a;
57c4762a1bSJed Brown   PetscScalar       *r;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A,&a));
629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(R,&r));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   r[0] = a[0] + (Omega*Omega)*u[0];
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A,&a));
689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(R,&r));
699566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(R));
709566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(R));
71c4762a1bSJed Brown   PetscFunctionReturn(0);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown PetscErrorCode Tangent1(TS ts,PetscReal t,Vec U,Vec A,PetscReal shiftA,Mat J,Mat P,void *ctx)
75c4762a1bSJed Brown {
76c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
77c4762a1bSJed Brown   PetscReal      Omega = user->Omega;
78c4762a1bSJed Brown   PetscReal      T = 0;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   PetscFunctionBegin;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   T = shiftA + (Omega*Omega);
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P,0,0,T,INSERT_VALUES));
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd  (P,MAT_FINAL_ASSEMBLY));
87c4762a1bSJed Brown   if (J != P) {
889566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
90c4762a1bSJed Brown   }
91c4762a1bSJed Brown   PetscFunctionReturn(0);
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown PetscErrorCode Residual2(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec R,void *ctx)
95c4762a1bSJed Brown {
96c4762a1bSJed Brown   UserParams         *user = (UserParams*)ctx;
97c4762a1bSJed Brown   PetscReal          Omega = user->Omega, Xi = user->Xi;
98c4762a1bSJed Brown   const PetscScalar *u,*v,*a;
99c4762a1bSJed Brown   PetscScalar       *r;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscFunctionBegin;
1029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V,&v));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A,&a));
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(R,&r));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   r[0] = a[0] + (2*Xi*Omega)*v[0] + (Omega*Omega)*u[0];
108c4762a1bSJed Brown 
1099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V,&v));
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A,&a));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(R,&r));
1139566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(R));
1149566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(R));
115c4762a1bSJed Brown   PetscFunctionReturn(0);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown PetscErrorCode Tangent2(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P,void *ctx)
119c4762a1bSJed Brown {
120c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
121c4762a1bSJed Brown   PetscReal      Omega = user->Omega, Xi = user->Xi;
122c4762a1bSJed Brown   PetscReal      T = 0;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   PetscFunctionBegin;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   T = shiftA + shiftV * (2*Xi*Omega) + (Omega*Omega);
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P,0,0,T,INSERT_VALUES));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd  (P,MAT_FINAL_ASSEMBLY));
131c4762a1bSJed Brown   if (J != P) {
1329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1339566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown int main(int argc, char *argv[])
139c4762a1bSJed Brown {
140c4762a1bSJed Brown   PetscMPIInt    size;
141c4762a1bSJed Brown   TS             ts;
142c4762a1bSJed Brown   Vec            R;
143c4762a1bSJed Brown   Mat            J;
144c4762a1bSJed Brown   Vec            U,V;
145c4762a1bSJed Brown   PetscScalar    *u,*v;
146c4762a1bSJed Brown   UserParams     user = {/*Omega=*/ 1, /*Xi=*/ 0, /*u0=*/ 1, /*,v0=*/ 0};
147c4762a1bSJed Brown 
148*327415f7SBarry Smith   PetscFunctionBeginUser;
1499566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
1509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1513c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
152c4762a1bSJed Brown 
153d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF,"","ex43 options","");
1549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-frequency","Natual frequency",__FILE__,user.Omega,&user.Omega,NULL));
1559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-damping","Damping coefficient",__FILE__,user.Xi,&user.Xi,NULL));
1569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-initial_u","Initial displacement",__FILE__,user.u0,&user.u0,NULL));
1579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-initial_v","Initial velocity",__FILE__,user.v0,&user.v0,NULL));
158d0609cedSBarry Smith   PetscOptionsEnd();
159c4762a1bSJed Brown 
1609566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
1619566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSALPHA2));
1629566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,5*(2*PETSC_PI)));
1639566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1649566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,0.01));
165c4762a1bSJed Brown 
1669566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,1,&R));
1679566063dSJacob Faibussowitsch   PetscCall(VecSetUp(R));
1689566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,1,1,NULL,&J));
1699566063dSJacob Faibussowitsch   PetscCall(MatSetUp(J));
170c4762a1bSJed Brown   if (user.Xi) {
1719566063dSJacob Faibussowitsch     PetscCall(TSSetI2Function(ts,R,Residual2,&user));
1729566063dSJacob Faibussowitsch     PetscCall(TSSetI2Jacobian(ts,J,J,Tangent2,&user));
173c4762a1bSJed Brown   } else {
1749566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts,R,Residual1,&user));
1759566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts,J,J,Tangent1,&user));
176c4762a1bSJed Brown   }
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
1789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1799566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts,Solution,&user));
180c4762a1bSJed Brown 
1819566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,1,&U));
1829566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,1,&V));
1839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(U,&u));
1849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(V,&v));
185c4762a1bSJed Brown   u[0] = user.u0;
186c4762a1bSJed Brown   v[0] = user.v0;
1879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(U,&u));
1889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(V,&v));
189c4762a1bSJed Brown 
1909566063dSJacob Faibussowitsch   PetscCall(TS2SetSolution(ts,U,V));
1919566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1929566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,NULL));
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V));
1969566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1979566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
198b122ec5aSJacob Faibussowitsch   return 0;
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown /*TEST
202c4762a1bSJed Brown 
203c4762a1bSJed Brown     test:
204c4762a1bSJed Brown       suffix: a
205c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_view
206c4762a1bSJed Brown       requires: !single
207c4762a1bSJed Brown 
208c4762a1bSJed Brown     test:
209c4762a1bSJed Brown       suffix: b
210c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_rtol 0 -ts_atol 1e-5 -ts_adapt_type basic -ts_adapt_monitor
211c4762a1bSJed Brown       requires: !single
212c4762a1bSJed Brown 
213c4762a1bSJed Brown TEST*/
214