xref: /petsc/src/ts/tutorials/ex43.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Single-DOF oscillator formulated as a second-order system.\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscts.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown typedef struct {
6*c4762a1bSJed Brown   PetscReal Omega;   /* natural frequency */
7*c4762a1bSJed Brown   PetscReal Xi;      /* damping coefficient  */
8*c4762a1bSJed Brown   PetscReal u0,v0;   /* initial conditions */
9*c4762a1bSJed Brown } UserParams;
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown static void Exact(PetscReal t,
12*c4762a1bSJed Brown                   PetscReal omega,PetscReal xi,PetscReal u0,PetscReal v0,
13*c4762a1bSJed Brown                   PetscReal *ut,PetscReal *vt)
14*c4762a1bSJed Brown {
15*c4762a1bSJed Brown #define Sin  PetscSinReal
16*c4762a1bSJed Brown #define Cos  PetscCosReal
17*c4762a1bSJed Brown #define Exp  PetscExpReal
18*c4762a1bSJed Brown #define Sqrt PetscSqrtReal
19*c4762a1bSJed Brown   PetscReal u,v;
20*c4762a1bSJed Brown   if (xi < 1) {
21*c4762a1bSJed Brown     PetscReal a  = xi*omega;
22*c4762a1bSJed Brown     PetscReal w  = Sqrt(1-xi*xi)*omega;
23*c4762a1bSJed Brown     PetscReal C1 = (v0 + a*u0)/w;
24*c4762a1bSJed Brown     PetscReal C2 = u0;
25*c4762a1bSJed Brown     u = Exp(-a*t) * (C1*Sin(w*t) + C2*Cos(w*t));
26*c4762a1bSJed Brown     v = (- a * Exp(-a*t) * (C1*Sin(w*t) + C2*Cos(w*t))
27*c4762a1bSJed Brown          + w * Exp(-a*t) * (C1*Cos(w*t) - C2*Sin(w*t)));
28*c4762a1bSJed Brown   } else if (xi > 1) {
29*c4762a1bSJed Brown     PetscReal w  = Sqrt(xi*xi-1)*omega;
30*c4762a1bSJed Brown     PetscReal C1 = (w*u0 + xi*u0 + v0)/(2*w);
31*c4762a1bSJed Brown     PetscReal C2 = (w*u0 - xi*u0 - v0)/(2*w);
32*c4762a1bSJed Brown     u = C1*Exp((-xi+w)*t) + C2*Exp((-xi-w)*t);
33*c4762a1bSJed Brown     v = C1*(-xi+w)*Exp((-xi+w)*t) + C2*(-xi-w)*Exp((-xi-w)*t);
34*c4762a1bSJed Brown   } else {
35*c4762a1bSJed Brown     PetscReal a  = xi*omega;
36*c4762a1bSJed Brown     PetscReal C1 = v0 + a*u0;
37*c4762a1bSJed Brown     PetscReal C2 = u0;
38*c4762a1bSJed Brown     u = (C1*t + C2) * Exp(-a*t);
39*c4762a1bSJed Brown     v = (C1 - a*(C1*t + C2)) * Exp(-a*t);
40*c4762a1bSJed Brown   }
41*c4762a1bSJed Brown   if (ut) *ut = u;
42*c4762a1bSJed Brown   if (vt) *vt = v;
43*c4762a1bSJed Brown }
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown PetscErrorCode Solution(TS ts,PetscReal t,Vec X,void *ctx)
46*c4762a1bSJed Brown {
47*c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
48*c4762a1bSJed Brown   PetscReal      u,v;
49*c4762a1bSJed Brown   PetscScalar    *x;
50*c4762a1bSJed Brown   PetscErrorCode ierr;
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   PetscFunctionBegin;
53*c4762a1bSJed Brown   Exact(t,user->Omega,user->Xi,user->u0,user->v0,&u,&v);
54*c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
55*c4762a1bSJed Brown   x[0] = (PetscScalar)u;
56*c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
57*c4762a1bSJed Brown   PetscFunctionReturn(0);
58*c4762a1bSJed Brown }
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown PetscErrorCode Residual1(TS ts,PetscReal t,Vec U,Vec A,Vec R,void *ctx)
61*c4762a1bSJed Brown {
62*c4762a1bSJed Brown   UserParams        *user = (UserParams*)ctx;
63*c4762a1bSJed Brown   PetscReal         Omega = user->Omega;
64*c4762a1bSJed Brown   const PetscScalar *u,*a;
65*c4762a1bSJed Brown   PetscScalar       *r;
66*c4762a1bSJed Brown   PetscErrorCode    ierr;
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   PetscFunctionBegin;
69*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = VecGetArrayRead(A,&a);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   r[0] = a[0] + (Omega*Omega)*u[0];
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(A,&a);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = VecAssemblyBegin(R);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = VecAssemblyEnd(R);CHKERRQ(ierr);
80*c4762a1bSJed Brown   PetscFunctionReturn(0);
81*c4762a1bSJed Brown }
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown PetscErrorCode Tangent1(TS ts,PetscReal t,Vec U,Vec A,PetscReal shiftA,Mat J,Mat P,void *ctx)
84*c4762a1bSJed Brown {
85*c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
86*c4762a1bSJed Brown   PetscReal      Omega = user->Omega;
87*c4762a1bSJed Brown   PetscReal      T = 0;
88*c4762a1bSJed Brown   PetscErrorCode ierr;
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   PetscFunctionBegin;
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   T = shiftA + (Omega*Omega);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   ierr = MatSetValue(P,0,0,T,INSERT_VALUES);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
97*c4762a1bSJed Brown   if (J != P) {
98*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
100*c4762a1bSJed Brown   }
101*c4762a1bSJed Brown   PetscFunctionReturn(0);
102*c4762a1bSJed Brown }
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown PetscErrorCode Residual2(TS ts,PetscReal t,Vec U,Vec V,Vec A,Vec R,void *ctx)
105*c4762a1bSJed Brown {
106*c4762a1bSJed Brown   UserParams         *user = (UserParams*)ctx;
107*c4762a1bSJed Brown   PetscReal          Omega = user->Omega, Xi = user->Xi;
108*c4762a1bSJed Brown   const PetscScalar *u,*v,*a;
109*c4762a1bSJed Brown   PetscScalar       *r;
110*c4762a1bSJed Brown   PetscErrorCode    ierr;
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown   PetscFunctionBegin;
113*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = VecGetArrayRead(V,&v);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = VecGetArrayRead(A,&a);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown   r[0] = a[0] + (2*Xi*Omega)*v[0] + (Omega*Omega)*u[0];
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(V,&v);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(A,&a);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = VecAssemblyBegin(R);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = VecAssemblyEnd  (R);CHKERRQ(ierr);
126*c4762a1bSJed Brown   PetscFunctionReturn(0);
127*c4762a1bSJed Brown }
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown PetscErrorCode Tangent2(TS ts,PetscReal t,Vec U,Vec V,Vec A,PetscReal shiftV,PetscReal shiftA,Mat J,Mat P,void *ctx)
130*c4762a1bSJed Brown {
131*c4762a1bSJed Brown   UserParams     *user = (UserParams*)ctx;
132*c4762a1bSJed Brown   PetscReal      Omega = user->Omega, Xi = user->Xi;
133*c4762a1bSJed Brown   PetscReal      T = 0;
134*c4762a1bSJed Brown   PetscErrorCode ierr;
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   PetscFunctionBegin;
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   T = shiftA + shiftV * (2*Xi*Omega) + (Omega*Omega);
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown   ierr = MatSetValue(P,0,0,T,INSERT_VALUES);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143*c4762a1bSJed Brown   if (J != P) {
144*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146*c4762a1bSJed Brown   }
147*c4762a1bSJed Brown   PetscFunctionReturn(0);
148*c4762a1bSJed Brown }
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown int main(int argc, char *argv[])
151*c4762a1bSJed Brown {
152*c4762a1bSJed Brown   PetscMPIInt    size;
153*c4762a1bSJed Brown   TS             ts;
154*c4762a1bSJed Brown   Vec            R;
155*c4762a1bSJed Brown   Mat            J;
156*c4762a1bSJed Brown   Vec            U,V;
157*c4762a1bSJed Brown   PetscScalar    *u,*v;
158*c4762a1bSJed Brown   UserParams     user = {/*Omega=*/ 1, /*Xi=*/ 0, /*u0=*/ 1, /*,v0=*/ 0};
159*c4762a1bSJed Brown   PetscErrorCode ierr;
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
162*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
163*c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_SELF,"","ex43 options","");CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = PetscOptionsReal("-frequency","Natual frequency",__FILE__,user.Omega,&user.Omega,NULL);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = PetscOptionsReal("-damping","Damping coefficient",__FILE__,user.Xi,&user.Xi,NULL);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = PetscOptionsReal("-initial_u","Initial displacement",__FILE__,user.u0,&user.u0,NULL);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = PetscOptionsReal("-initial_v","Initial velocity",__FILE__,user.v0,&user.v0,NULL);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = TSSetType(ts,TSALPHA2);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,5*(2*PETSC_PI));CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr);
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&R);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = VecSetUp(R);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = MatCreateSeqDense(PETSC_COMM_SELF,1,1,NULL,&J);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = MatSetUp(J);CHKERRQ(ierr);
182*c4762a1bSJed Brown   if (user.Xi) {
183*c4762a1bSJed Brown     ierr = TSSetI2Function(ts,R,Residual2,&user);CHKERRQ(ierr);
184*c4762a1bSJed Brown     ierr = TSSetI2Jacobian(ts,J,J,Tangent2,&user);CHKERRQ(ierr);
185*c4762a1bSJed Brown   } else {
186*c4762a1bSJed Brown     ierr = TSSetIFunction(ts,R,Residual1,&user);CHKERRQ(ierr);
187*c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,J,J,Tangent1,&user);CHKERRQ(ierr);
188*c4762a1bSJed Brown   }
189*c4762a1bSJed Brown   ierr = VecDestroy(&R);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = TSSetSolutionFunction(ts,Solution,&user);CHKERRQ(ierr);
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&U);CHKERRQ(ierr);
194*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_SELF,1,&V);CHKERRQ(ierr);
195*c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
196*c4762a1bSJed Brown   ierr = VecGetArray(V,&v);CHKERRQ(ierr);
197*c4762a1bSJed Brown   u[0] = user.u0;
198*c4762a1bSJed Brown   v[0] = user.v0;
199*c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = VecRestoreArray(V,&v);CHKERRQ(ierr);
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown   ierr = TS2SetSolution(ts,U,V);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
204*c4762a1bSJed Brown   ierr = TSSolve(ts,NULL);CHKERRQ(ierr);
205*c4762a1bSJed Brown 
206*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
207*c4762a1bSJed Brown   ierr = VecDestroy(&V);CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
209*c4762a1bSJed Brown   ierr = PetscFinalize();
210*c4762a1bSJed Brown   return ierr;
211*c4762a1bSJed Brown }
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown /*TEST
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown     test:
216*c4762a1bSJed Brown       suffix: a
217*c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_view
218*c4762a1bSJed Brown       requires: !single
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown     test:
221*c4762a1bSJed Brown       suffix: b
222*c4762a1bSJed Brown       args: -ts_max_steps 10 -ts_rtol 0 -ts_atol 1e-5 -ts_adapt_type basic -ts_adapt_monitor
223*c4762a1bSJed Brown       requires: !single
224*c4762a1bSJed Brown 
225*c4762a1bSJed Brown TEST*/
226