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