xref: /petsc/src/ts/tutorials/ex43.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
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   PetscErrorCode ierr;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   PetscFunctionBegin;
46c4762a1bSJed Brown   Exact(t,user->Omega,user->Xi,user->u0,user->v0,&u,&v);
47c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
48c4762a1bSJed Brown   x[0] = (PetscScalar)u;
49c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
50c4762a1bSJed Brown   PetscFunctionReturn(0);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown PetscErrorCode Residual1(TS ts,PetscReal t,Vec U,Vec A,Vec R,void *ctx)
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   UserParams        *user = (UserParams*)ctx;
56c4762a1bSJed Brown   PetscReal         Omega = user->Omega;
57c4762a1bSJed Brown   const PetscScalar *u,*a;
58c4762a1bSJed Brown   PetscScalar       *r;
59c4762a1bSJed Brown   PetscErrorCode    ierr;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   PetscFunctionBegin;
62c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = VecGetArrayRead(A,&a);CHKERRQ(ierr);
64303a5415SBarry Smith   ierr = VecGetArrayWrite(R,&r);CHKERRQ(ierr);
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   r[0] = a[0] + (Omega*Omega)*u[0];
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = VecRestoreArrayRead(A,&a);CHKERRQ(ierr);
70303a5415SBarry Smith   ierr = VecRestoreArrayWrite(R,&r);CHKERRQ(ierr);
71c4762a1bSJed Brown   ierr = VecAssemblyBegin(R);CHKERRQ(ierr);
72c4762a1bSJed Brown   ierr = VecAssemblyEnd(R);CHKERRQ(ierr);
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown PetscErrorCode Tangent1(TS ts,PetscReal t,Vec U,Vec A,PetscReal shiftA,Mat J,Mat P,void *ctx)
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
79c4762a1bSJed Brown   PetscReal      Omega = user->Omega;
80c4762a1bSJed Brown   PetscReal      T = 0;
81c4762a1bSJed Brown   PetscErrorCode ierr;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   PetscFunctionBegin;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   T = shiftA + (Omega*Omega);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   ierr = MatSetValue(P,0,0,T,INSERT_VALUES);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = MatAssemblyEnd  (P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90c4762a1bSJed Brown   if (J != P) {
91c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown PetscErrorCode Residual2(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec R,void *ctx)
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   UserParams         *user = (UserParams*)ctx;
100c4762a1bSJed Brown   PetscReal          Omega = user->Omega, Xi = user->Xi;
101c4762a1bSJed Brown   const PetscScalar *u,*v,*a;
102c4762a1bSJed Brown   PetscScalar       *r;
103c4762a1bSJed Brown   PetscErrorCode    ierr;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   PetscFunctionBegin;
106c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = VecGetArrayRead(V,&v);CHKERRQ(ierr);
108c4762a1bSJed Brown   ierr = VecGetArrayRead(A,&a);CHKERRQ(ierr);
109303a5415SBarry Smith   ierr = VecGetArrayWrite(R,&r);CHKERRQ(ierr);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   r[0] = a[0] + (2*Xi*Omega)*v[0] + (Omega*Omega)*u[0];
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = VecRestoreArrayRead(V,&v);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = VecRestoreArrayRead(A,&a);CHKERRQ(ierr);
116303a5415SBarry Smith   ierr = VecRestoreArrayWrite(R,&r);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = VecAssemblyBegin(R);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = VecAssemblyEnd(R);CHKERRQ(ierr);
119c4762a1bSJed Brown   PetscFunctionReturn(0);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown PetscErrorCode Tangent2(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P,void *ctx)
123c4762a1bSJed Brown {
124c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
125c4762a1bSJed Brown   PetscReal      Omega = user->Omega, Xi = user->Xi;
126c4762a1bSJed Brown   PetscReal      T = 0;
127c4762a1bSJed Brown   PetscErrorCode ierr;
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   PetscFunctionBegin;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   T = shiftA + shiftV * (2*Xi*Omega) + (Omega*Omega);
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   ierr = MatSetValue(P,0,0,T,INSERT_VALUES);CHKERRQ(ierr);
134c4762a1bSJed Brown   ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135c4762a1bSJed Brown   ierr = MatAssemblyEnd  (P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136c4762a1bSJed Brown   if (J != P) {
137c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown   PetscFunctionReturn(0);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown int main(int argc, char *argv[])
144c4762a1bSJed Brown {
145c4762a1bSJed Brown   PetscMPIInt    size;
146c4762a1bSJed Brown   TS             ts;
147c4762a1bSJed Brown   Vec            R;
148c4762a1bSJed Brown   Mat            J;
149c4762a1bSJed Brown   Vec            U,V;
150c4762a1bSJed Brown   PetscScalar    *u,*v;
151c4762a1bSJed Brown   UserParams     user = {/*Omega=*/ 1, /*Xi=*/ 0, /*u0=*/ 1, /*,v0=*/ 0};
152c4762a1bSJed Brown   PetscErrorCode ierr;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
155*ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
156c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_SELF,"","ex43 options","");CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = PetscOptionsReal("-frequency","Natual frequency",__FILE__,user.Omega,&user.Omega,NULL);CHKERRQ(ierr);
160c4762a1bSJed Brown   ierr = PetscOptionsReal("-damping","Damping coefficient",__FILE__,user.Xi,&user.Xi,NULL);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = PetscOptionsReal("-initial_u","Initial displacement",__FILE__,user.u0,&user.u0,NULL);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = PetscOptionsReal("-initial_v","Initial velocity",__FILE__,user.v0,&user.v0,NULL);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr);
166c4762a1bSJed Brown   ierr = TSSetType(ts,TSALPHA2);CHKERRQ(ierr);
167c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,5*(2*PETSC_PI));CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr);
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&R);CHKERRQ(ierr);
172c4762a1bSJed Brown   ierr = VecSetUp(R);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = MatCreateSeqDense(PETSC_COMM_SELF,1,1,NULL,&J);CHKERRQ(ierr);
174c4762a1bSJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
175c4762a1bSJed Brown   if (user.Xi) {
176c4762a1bSJed Brown     ierr = TSSetI2Function(ts,R,Residual2,&user);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = TSSetI2Jacobian(ts,J,J,Tangent2,&user);CHKERRQ(ierr);
178c4762a1bSJed Brown   } else {
179c4762a1bSJed Brown     ierr = TSSetIFunction(ts,R,Residual1,&user);CHKERRQ(ierr);
180c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,J,J,Tangent1,&user);CHKERRQ(ierr);
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown   ierr = VecDestroy(&R);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
184c4762a1bSJed Brown   ierr = TSSetSolutionFunction(ts,Solution,&user);CHKERRQ(ierr);
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&U);CHKERRQ(ierr);
187c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&V);CHKERRQ(ierr);
188303a5415SBarry Smith   ierr = VecGetArrayWrite(U,&u);CHKERRQ(ierr);
189303a5415SBarry Smith   ierr = VecGetArrayWrite(V,&v);CHKERRQ(ierr);
190c4762a1bSJed Brown   u[0] = user.u0;
191c4762a1bSJed Brown   v[0] = user.v0;
192303a5415SBarry Smith   ierr = VecRestoreArrayWrite(U,&u);CHKERRQ(ierr);
193303a5415SBarry Smith   ierr = VecRestoreArrayWrite(V,&v);CHKERRQ(ierr);
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   ierr = TS2SetSolution(ts,U,V);CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = TSSolve(ts,NULL);CHKERRQ(ierr);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
200c4762a1bSJed Brown   ierr = VecDestroy(&V);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = PetscFinalize();
203c4762a1bSJed Brown   return ierr;
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*TEST
207c4762a1bSJed Brown 
208c4762a1bSJed Brown     test:
209c4762a1bSJed Brown       suffix: a
210c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_view
211c4762a1bSJed Brown       requires: !single
212c4762a1bSJed Brown 
213c4762a1bSJed Brown     test:
214c4762a1bSJed Brown       suffix: b
215c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_rtol 0 -ts_atol 1e-5 -ts_adapt_type basic -ts_adapt_monitor
216c4762a1bSJed Brown       requires: !single
217c4762a1bSJed Brown 
218c4762a1bSJed Brown TEST*/
219