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