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) = -6.6e-01,
17    and
18        mu = 10^6.
19    This is a nonlinear equation.
20 
21    This is a copy and modification of ex20.c to exactly match a test
22    problem that comes with the Radau5 integrator package.
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 
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void * ctx)35 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
36 {
37   PetscErrorCode    ierr;
38   User              user = (User)ctx;
39   const PetscScalar *x,*xdot;
40   PetscScalar       *f;
41 
42   PetscFunctionBeginUser;
43   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
44   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
45   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
46   f[0] = xdot[0] - x[1];
47   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
48   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
49   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
50   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
51   PetscFunctionReturn(0);
52 }
53 
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void * ctx)54 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
55 {
56   PetscErrorCode    ierr;
57   User              user     = (User)ctx;
58   PetscInt          rowcol[] = {0,1};
59   const PetscScalar *x;
60   PetscScalar       J[2][2];
61 
62   PetscFunctionBeginUser;
63   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
64   J[0][0] = a;     J[0][1] = -1.0;
65   J[1][0] = user->mu*(1.0 + 2.0*x[0]*x[1]);   J[1][1] = a - user->mu*(1.0-x[0]*x[0]);
66   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
67   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
68 
69   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71   if (A != B) {
72     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
73     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74   }
75   PetscFunctionReturn(0);
76 }
77 
78 
main(int argc,char ** argv)79 int main(int argc,char **argv)
80 {
81   TS             ts;            /* nonlinear solver */
82   Vec            x;             /* solution, residual vectors */
83   Mat            A;             /* Jacobian matrix */
84   PetscInt       steps;
85   PetscReal      ftime   = 2;
86   PetscScalar    *x_ptr;
87   PetscMPIInt    size;
88   struct _n_User user;
89   PetscErrorCode ierr;
90 
91   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92      Initialize program
93      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
95   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
96   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
97 
98   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99     Set runtime options
100     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101   user.next_output = 0.0;
102   user.mu          = 1.0e6;
103   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
104   ierr = PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);CHKERRQ(ierr);
105   ierr = PetscOptionsEnd();CHKERRQ(ierr);
106 
107   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108     Create necessary matrix and vectors, solve same ODE on every process
109     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
111   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
112   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
113   ierr = MatSetUp(A);CHKERRQ(ierr);
114 
115   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
116 
117   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118      Create timestepping solver context
119      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
121   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
122   ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
123   ierr = TSSetIJacobian(ts,A,A,IJacobian,&user);CHKERRQ(ierr);
124 
125   ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr);
126   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
127   ierr = TSSetTolerances(ts,1.e-4,NULL,1.e-4,NULL);CHKERRQ(ierr);
128   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129      Set initial conditions
130    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131   ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr);
132   x_ptr[0] = 2.0;   x_ptr[1] = -6.6e-01;
133   ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr);
134   ierr = TSSetTimeStep(ts,.000001);CHKERRQ(ierr);
135 
136   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137      Set runtime options
138    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
140 
141   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142      Solve nonlinear system
143      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144   ierr = TSSolve(ts,x);CHKERRQ(ierr);
145   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
146   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
147   ierr = PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);CHKERRQ(ierr);
148   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
149 
150   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151      Free work space.  All PETSc objects should be destroyed when they
152      are no longer needed.
153    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154   ierr = MatDestroy(&A);CHKERRQ(ierr);
155   ierr = VecDestroy(&x);CHKERRQ(ierr);
156   ierr = TSDestroy(&ts);CHKERRQ(ierr);
157 
158   ierr = PetscFinalize();
159   return(ierr);
160 }
161 
162 /*TEST
163 
164     build:
165       requires: double !complex !define(PETSC_USE_64BIT_INDICES) radau5
166 
167     test:
168       args: -ts_monitor_solution -ts_type radau5
169 
170 TEST*/
171