xref: /petsc/src/ts/tutorials/ex16.c (revision efbe7e8a80d07327753dbe0b33efee01e046af3f)
1 
2 static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\
3 Input parameters include:\n\
4       -mu : stiffness parameter\n\n";
5 
6 /*
7    Concepts: TS^time-dependent nonlinear problems
8    Concepts: TS^van der Pol equation
9    Processors: 1
10 */
11 /* ------------------------------------------------------------------------
12 
13    This program solves the van der Pol equation
14        y'' - \mu ((1-y^2)*y' - y) = 0        (1)
15    on the domain 0 <= x <= 1, with the boundary conditions
16        y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
17    This is a nonlinear equation. The well prepared initial condition gives errors that are not dominated by the first few steps of the method when \mu is large.
18 
19    Notes:
20    This code demonstrates the TS solver interface to two variants of
21    linear problems, u_t = f(u,t), namely turning (1) into a system of
22    first order differential equations,
23 
24    [ y' ] = [          z            ]
25    [ z' ]   [ \mu ((1 - y^2) z - y) ]
26 
27    which then we can write as a vector equation
28 
29    [ u_1' ] = [             u_2           ]  (2)
30    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]
31 
32    which is now in the desired form of u_t = f(u,t). One way that we
33    can split f(u,t) in (2) is to split by component,
34 
35    [ u_1' ] = [ u_2 ] + [            0                ]
36    [ u_2' ]   [  0  ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]
37 
38    where
39 
40    [ G(u,t) ] = [ u_2 ]
41                 [  0  ]
42 
43    and
44 
45    [ F(u',u,t) ] = [ u_1' ] - [            0                ]
46                    [ u_2' ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]
47 
48    Using the definition of the Jacobian of F (from the PETSc user manual),
49    in the equation F(u',u,t) = G(u,t),
50 
51               dF   dF
52    J(F) = a * -- - --
53               du'  du
54 
55    where d is the partial derivative. In this example,
56 
57    dF   [ 1 ; 0 ]
58    -- = [       ]
59    du'  [ 0 ; 1 ]
60 
61    dF   [       0             ;         0        ]
62    -- = [                                        ]
63    du   [ -\mu (2*u_1*u_2 + 1);  \mu (1 - u_1^2) ]
64 
65    Hence,
66 
67           [      a             ;          0          ]
68    J(F) = [                                          ]
69           [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ]
70 
71   ------------------------------------------------------------------------- */
72 
73 #include <petscts.h>
74 
75 typedef struct _n_User *User;
76 struct _n_User {
77   PetscReal mu;
78   PetscBool imex;
79   PetscReal next_output;
80 };
81 
82 /*
83 *  User-defined routines
84 */
85 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
86 {
87   PetscErrorCode    ierr;
88   User              user = (User)ctx;
89   PetscScalar       *f;
90   const PetscScalar *x;
91 
92   PetscFunctionBeginUser;
93   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
94   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
95   f[0] = (user->imex ? x[1] : 0);
96   f[1] = 0.0;
97   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
98   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
99   PetscFunctionReturn(0);
100 }
101 
102 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
103 {
104   PetscErrorCode    ierr;
105   User              user = (User)ctx;
106   const PetscScalar *x,*xdot;
107   PetscScalar       *f;
108 
109   PetscFunctionBeginUser;
110   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
111   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
112   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
113   f[0] = xdot[0] + (user->imex ? 0 : x[1]);
114   f[1] = xdot[1] - user->mu*((1. - x[0]*x[0])*x[1] - x[0]);
115   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
116   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
117   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
118   PetscFunctionReturn(0);
119 }
120 
121 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
122 {
123   PetscErrorCode    ierr;
124   User              user     = (User)ctx;
125   PetscReal         mu       = user->mu;
126   PetscInt          rowcol[] = {0,1};
127   const PetscScalar *x;
128   PetscScalar       J[2][2];
129 
130   PetscFunctionBeginUser;
131   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
132   J[0][0] = a;                    J[0][1] = (user->imex ? 0 : 1.);
133   J[1][0] = mu*(2.*x[0]*x[1]+1.);   J[1][1] = a - mu*(1. - x[0]*x[0]);
134   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
135   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
136 
137   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139   if (A != B) {
140     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
141     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
142   }
143   PetscFunctionReturn(0);
144 }
145 
146 static PetscErrorCode RegisterMyARK2(void)
147 {
148   PetscErrorCode ierr;
149 
150   PetscFunctionBeginUser;
151   {
152     const PetscReal
153       A[3][3] = {{0,0,0},
154                  {0.41421356237309504880,0,0},
155                  {0.75,0.25,0}},
156       At[3][3] = {{0,0,0},
157                   {0.12132034355964257320,0.29289321881345247560,0},
158                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}},
159       *bembedt = NULL,*bembed = NULL;
160     ierr = TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,bembedt,bembed,0,NULL,NULL);CHKERRQ(ierr);
161   }
162   PetscFunctionReturn(0);
163 }
164 
165 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
166 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
167 {
168   PetscErrorCode    ierr;
169   const PetscScalar *x;
170   PetscReal         tfinal, dt;
171   User              user = (User)ctx;
172   Vec               interpolatedX;
173 
174   PetscFunctionBeginUser;
175   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
176   ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr);
177 
178   while (user->next_output <= t && user->next_output <= tfinal) {
179     ierr = VecDuplicate(X,&interpolatedX);CHKERRQ(ierr);
180     ierr = TSInterpolate(ts,user->next_output,interpolatedX);CHKERRQ(ierr);
181     ierr = VecGetArrayRead(interpolatedX,&x);CHKERRQ(ierr);
182     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",user->next_output,step,t,dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]));CHKERRQ(ierr);
183     ierr = VecRestoreArrayRead(interpolatedX,&x);CHKERRQ(ierr);
184     ierr = VecDestroy(&interpolatedX);CHKERRQ(ierr);
185 
186     user->next_output += 0.1;
187   }
188   PetscFunctionReturn(0);
189 }
190 
191 int main(int argc,char **argv)
192 {
193   TS             ts;            /* nonlinear solver */
194   Vec            x;             /* solution, residual vectors */
195   Mat            A;             /* Jacobian matrix */
196   PetscInt       steps;
197   PetscReal      ftime = 0.5;
198   PetscBool      monitor = PETSC_FALSE;
199   PetscScalar    *x_ptr;
200   PetscMPIInt    size;
201   struct _n_User user;
202   PetscErrorCode ierr;
203 
204   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205      Initialize program
206      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
208   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
209   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
210 
211   ierr = RegisterMyARK2();CHKERRQ(ierr);
212 
213   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214     Set runtime options
215     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216   user.mu          = 1000.0;
217   user.imex        = PETSC_TRUE;
218   user.next_output = 0.0;
219 
220   ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr);
221   ierr = PetscOptionsGetBool(NULL,NULL,"-imex",&user.imex,NULL);CHKERRQ(ierr);
222   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
223 
224   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225     Create necessary matrix and vectors, solve same ODE on every process
226     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
228   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
229   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
230   ierr = MatSetUp(A);CHKERRQ(ierr);
231   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
232 
233   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234      Create timestepping solver context
235      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
237   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
238   ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
239   ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
240   ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr);
241   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
242   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
243   if (monitor) {
244     ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr);
245   }
246 
247   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248      Set initial conditions
249    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250   ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr);
251   x_ptr[0] = 2.0;
252   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
253   ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr);
254   ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr);
255 
256   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257      Set runtime options
258    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
259   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
260 
261   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262      Solve nonlinear system
263      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264   ierr = TSSolve(ts,x);CHKERRQ(ierr);
265   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
266   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
267   ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);CHKERRQ(ierr);
268   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
269 
270   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271      Free work space.  All PETSc objects should be destroyed when they
272      are no longer needed.
273    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274   ierr = MatDestroy(&A);CHKERRQ(ierr);
275   ierr = VecDestroy(&x);CHKERRQ(ierr);
276   ierr = TSDestroy(&ts);CHKERRQ(ierr);
277 
278   ierr = PetscFinalize();
279   return ierr;
280 }
281 
282 /*TEST
283 
284     test:
285       args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none
286       requires: !single
287 
288 TEST*/
289