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