xref: /petsc/src/ts/tutorials/ex49.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
7c4762a1bSJed Brown    This program solves the van der Pol DAE ODE equivalent
8c4762a1bSJed Brown        y' = z                 (1)
9c4762a1bSJed Brown        z' = mu[(1-y^2)z-y]
10c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
11c4762a1bSJed Brown        y(0) = 2, y'(0) = -6.6e-01,
12c4762a1bSJed Brown    and
13c4762a1bSJed Brown        mu = 10^6.
14c4762a1bSJed Brown    This is a nonlinear equation.
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    This is a copy and modification of ex20.c to exactly match a test
17c4762a1bSJed Brown    problem that comes with the Radau5 integrator package.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ------------------------------------------------------------------------- */
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown typedef struct _n_User *User;
24c4762a1bSJed Brown struct _n_User {
25c4762a1bSJed Brown   PetscReal mu;
26c4762a1bSJed Brown   PetscReal next_output;
27c4762a1bSJed Brown };
28c4762a1bSJed Brown 
29c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   User              user = (User)ctx;
32c4762a1bSJed Brown   const PetscScalar *x,*xdot;
33c4762a1bSJed Brown   PetscScalar       *f;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot,&xdot));
389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
39c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
40c4762a1bSJed Brown   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
44c4762a1bSJed Brown   PetscFunctionReturn(0);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
48c4762a1bSJed Brown {
49c4762a1bSJed Brown   User              user     = (User)ctx;
50c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
51c4762a1bSJed Brown   const PetscScalar *x;
52c4762a1bSJed Brown   PetscScalar       J[2][2];
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   PetscFunctionBeginUser;
559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
56c4762a1bSJed Brown   J[0][0] = a;     J[0][1] = -1.0;
57c4762a1bSJed Brown   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]);
589566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
63c4762a1bSJed Brown   if (A != B) {
649566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
659566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown   PetscFunctionReturn(0);
68c4762a1bSJed Brown }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown int main(int argc,char **argv)
71c4762a1bSJed Brown {
72c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
73c4762a1bSJed Brown   Vec            x;             /* solution, residual vectors */
74c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
75c4762a1bSJed Brown   PetscInt       steps;
76c4762a1bSJed Brown   PetscReal      ftime   = 2;
77c4762a1bSJed Brown   PetscScalar    *x_ptr;
78c4762a1bSJed Brown   PetscMPIInt    size;
79c4762a1bSJed Brown   struct _n_User user;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82c4762a1bSJed Brown      Initialize program
83c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
849566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
863c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4762a1bSJed Brown     Set runtime options
90c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91c4762a1bSJed Brown   user.next_output = 0.0;
92c4762a1bSJed Brown   user.mu          = 1.0e6;
93*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL));
95*d0609cedSBarry Smith   PetscOptionsEnd();
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
99c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1009566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
1019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
1029566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1039566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&x,NULL));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown      Create timestepping solver context
109c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1109566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1119566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBEULER));
1129566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,IFunction,&user));
1139566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,IJacobian,&user));
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,ftime));
1169566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1179566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(ts,1.e-4,NULL,1.e-4,NULL));
118c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown      Set initial conditions
120c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x,&x_ptr));
122c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = -6.6e-01;
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x,&x_ptr));
1249566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.000001));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Set runtime options
128c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1299566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132c4762a1bSJed Brown      Solve nonlinear system
133c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1349566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
1359566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
1369566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
1379566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime));
1389566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
142c4762a1bSJed Brown      are no longer needed.
143c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1469566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
149*d0609cedSBarry Smith   return(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /*TEST
153c4762a1bSJed Brown 
154c4762a1bSJed Brown     build:
155dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) radau5
156c4762a1bSJed Brown 
157c4762a1bSJed Brown     test:
158c4762a1bSJed Brown       args: -ts_monitor_solution -ts_type radau5
159c4762a1bSJed Brown 
160c4762a1bSJed Brown TEST*/
161