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