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