Actual source code: ex49.c


  2: static char help[] = "Solves the van der Pol equation.\n\
  3: Input parameters include:\n";

  5: /* ------------------------------------------------------------------------

  7:    This program solves the van der Pol DAE ODE equivalent
  8:        y' = z                 (1)
  9:        z' = mu[(1-y^2)z-y]
 10:    on the domain 0 <= x <= 1, with the boundary conditions
 11:        y(0) = 2, y'(0) = -6.6e-01,
 12:    and
 13:        mu = 10^6.
 14:    This is a nonlinear equation.

 16:    This is a copy and modification of ex20.c to exactly match a test
 17:    problem that comes with the Radau5 integrator package.

 19:   ------------------------------------------------------------------------- */

 21: #include <petscts.h>

 23: typedef struct _n_User *User;
 24: struct _n_User {
 25:   PetscReal mu;
 26:   PetscReal next_output;
 27: };

 29: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 30: {
 31:   User              user = (User)ctx;
 32:   const PetscScalar *x,*xdot;
 33:   PetscScalar       *f;

 36:   VecGetArrayRead(X,&x);
 37:   VecGetArrayRead(Xdot,&xdot);
 38:   VecGetArray(F,&f);
 39:   f[0] = xdot[0] - x[1];
 40:   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
 41:   VecRestoreArrayRead(X,&x);
 42:   VecRestoreArrayRead(Xdot,&xdot);
 43:   VecRestoreArray(F,&f);
 44:   return 0;
 45: }

 47: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 48: {
 49:   User              user     = (User)ctx;
 50:   PetscInt          rowcol[] = {0,1};
 51:   const PetscScalar *x;
 52:   PetscScalar       J[2][2];

 55:   VecGetArrayRead(X,&x);
 56:   J[0][0] = a;     J[0][1] = -1.0;
 57:   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]);
 58:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 59:   VecRestoreArrayRead(X,&x);

 61:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 62:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 63:   if (A != B) {
 64:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 65:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 66:   }
 67:   return 0;
 68: }

 70: int main(int argc,char **argv)
 71: {
 72:   TS             ts;            /* nonlinear solver */
 73:   Vec            x;             /* solution, residual vectors */
 74:   Mat            A;             /* Jacobian matrix */
 75:   PetscInt       steps;
 76:   PetscReal      ftime   = 2;
 77:   PetscScalar    *x_ptr;
 78:   PetscMPIInt    size;
 79:   struct _n_User user;

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:      Initialize program
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 85:   PetscInitialize(&argc,&argv,NULL,help);
 86:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 90:     Set runtime options
 91:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 92:   user.next_output = 0.0;
 93:   user.mu          = 1.0e6;
 94:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
 95:   PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);
 96:   PetscOptionsEnd();

 98:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 99:     Create necessary matrix and vectors, solve same ODE on every process
100:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101:   MatCreate(PETSC_COMM_WORLD,&A);
102:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
103:   MatSetFromOptions(A);
104:   MatSetUp(A);

106:   MatCreateVecs(A,&x,NULL);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:      Create timestepping solver context
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111:   TSCreate(PETSC_COMM_WORLD,&ts);
112:   TSSetType(ts,TSBEULER);
113:   TSSetIFunction(ts,NULL,IFunction,&user);
114:   TSSetIJacobian(ts,A,A,IJacobian,&user);

116:   TSSetMaxTime(ts,ftime);
117:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
118:   TSSetTolerances(ts,1.e-4,NULL,1.e-4,NULL);
119:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120:      Set initial conditions
121:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122:   VecGetArray(x,&x_ptr);
123:   x_ptr[0] = 2.0;   x_ptr[1] = -6.6e-01;
124:   VecRestoreArray(x,&x_ptr);
125:   TSSetTimeStep(ts,.000001);

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:      Set runtime options
129:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:   TSSetFromOptions(ts);

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Solve nonlinear system
134:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   TSSolve(ts,x);
136:   TSGetSolveTime(ts,&ftime);
137:   TSGetStepNumber(ts,&steps);
138:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
139:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:      Free work space.  All PETSc objects should be destroyed when they
143:      are no longer needed.
144:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   MatDestroy(&A);
146:   VecDestroy(&x);
147:   TSDestroy(&ts);

149:   PetscFinalize();
150:   return(ierr);
151: }

153: /*TEST

155:     build:
156:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) radau5

158:     test:
159:       args: -ts_monitor_solution -ts_type radau5

161: TEST*/