xref: /petsc/src/ts/tutorials/ex20.c (revision e607c864f3adee196c0614801d5954446e31f8d0)
1 
2 static char help[] = "Solves the van der Pol equation.\n\
3 Input parameters include:\n";
4 
5 /*
6    Concepts: TS^time-dependent nonlinear problems
7    Concepts: TS^van der Pol equation DAE equivalent
8    Processors: 1
9 */
10 /* ------------------------------------------------------------------------
11 
12    This program solves the van der Pol DAE ODE equivalent
13        y' = z                 (1)
14        z' = \mu ((1-y^2)z-y)
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    and
18        \mu = 10^6 ( y'(0) ~ -0.6666665432100101).
19    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.
20 
21    Notes:
22    This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form.
23 
24   ------------------------------------------------------------------------- */
25 
26 #include <petscts.h>
27 
28 typedef struct _n_User *User;
29 struct _n_User {
30   PetscReal mu;
31   PetscReal next_output;
32 };
33 
34 /*
35    User-defined routines
36 */
37 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
38 {
39   PetscErrorCode    ierr;
40   User              user = (User)ctx;
41   PetscScalar       *f;
42   const PetscScalar *x;
43 
44   PetscFunctionBeginUser;
45   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
46   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
47   f[0] = x[1];
48   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
49   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
50   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
54 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
55 {
56   PetscErrorCode    ierr;
57   User              user = (User)ctx;
58   const PetscScalar *x,*xdot;
59   PetscScalar       *f;
60 
61   PetscFunctionBeginUser;
62   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
63   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
64   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
65   f[0] = xdot[0] - x[1];
66   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
67   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
68   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
69   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
70   PetscFunctionReturn(0);
71 }
72 
73 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
74 {
75   PetscErrorCode    ierr;
76   User              user     = (User)ctx;
77   PetscInt          rowcol[] = {0,1};
78   const PetscScalar *x;
79   PetscScalar       J[2][2];
80 
81   PetscFunctionBeginUser;
82   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
83   J[0][0] = a;     J[0][1] = -1.0;
84   J[1][0] = user->mu*(2.0*x[0]*x[1] + 1.0);   J[1][1] = a - user->mu*(1.0-x[0]*x[0]);
85   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
86   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
87 
88   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90   if (A != B) {
91     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
98 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
99 {
100   PetscErrorCode    ierr;
101   const PetscScalar *x;
102   PetscReal         tfinal, dt;
103   User              user = (User)ctx;
104   Vec               interpolatedX;
105 
106   PetscFunctionBeginUser;
107   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
108   ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr);
109 
110   while (user->next_output <= t && user->next_output <= tfinal) {
111     ierr = VecDuplicate(X,&interpolatedX);CHKERRQ(ierr);
112     ierr = TSInterpolate(ts,user->next_output,interpolatedX);CHKERRQ(ierr);
113     ierr = VecGetArrayRead(interpolatedX,&x);CHKERRQ(ierr);
114     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
115                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
116                        (double)PetscRealPart(x[1]));CHKERRQ(ierr);
117     ierr = VecRestoreArrayRead(interpolatedX,&x);CHKERRQ(ierr);
118     ierr = VecDestroy(&interpolatedX);CHKERRQ(ierr);
119     user->next_output += 0.1;
120   }
121   PetscFunctionReturn(0);
122 }
123 
124 int main(int argc,char **argv)
125 {
126   TS             ts;            /* nonlinear solver */
127   Vec            x;             /* solution, residual vectors */
128   Mat            A;             /* Jacobian matrix */
129   PetscInt       steps;
130   PetscReal      ftime = 0.5;
131   PetscBool      monitor = PETSC_FALSE,implicitform = PETSC_TRUE;
132   PetscScalar    *x_ptr;
133   PetscMPIInt    size;
134   struct _n_User user;
135   PetscErrorCode ierr;
136 
137   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138      Initialize program
139      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
141   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
142   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
143 
144   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145     Set runtime options
146     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147   user.next_output = 0.0;
148   user.mu          = 1.0e3;
149   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
150   ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr);
151   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);CHKERRQ(ierr);
152   ierr = PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);CHKERRQ(ierr);
153   ierr = PetscOptionsEnd();CHKERRQ(ierr);
154 
155   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156     Create necessary matrix and vectors, solve same ODE on every process
157     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
159   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
160   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
161   ierr = MatSetUp(A);CHKERRQ(ierr);
162 
163   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
164 
165   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166      Create timestepping solver context
167      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
169   if (implicitform) {
170     ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
171     ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr);
172     ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
173   } else {
174     ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
175     ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
176   }
177   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
178   ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr);
179   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
180   if (monitor) {
181     ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr);
182   }
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185      Set initial conditions
186    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr);
188   x_ptr[0] = 2.0;
189   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
190   ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr);
191 
192   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193      Set runtime options
194    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
196 
197   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198      Solve nonlinear system
199      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200   ierr = TSSolve(ts,x);CHKERRQ(ierr);
201   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
202   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
203   ierr = PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);CHKERRQ(ierr);
204   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
205 
206   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207      Free work space.  All PETSc objects should be destroyed when they
208      are no longer needed.
209    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210   ierr = MatDestroy(&A);CHKERRQ(ierr);
211   ierr = VecDestroy(&x);CHKERRQ(ierr);
212   ierr = TSDestroy(&ts);CHKERRQ(ierr);
213 
214   ierr = PetscFinalize();
215   return(ierr);
216 }
217 
218 /*TEST
219 
220     test:
221       requires: !single
222       args: -mu 1e6
223 
224     test:
225       requires: !single
226       suffix: 2
227       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp
228 
229     test:
230       requires: !single
231       suffix: 3
232       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312
233 
234 TEST*/
235