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