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