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