Actual source code: ex49.c


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

  5: /*
  6:    Concepts: TS^time-dependent nonlinear problems
  7:    Concepts: TS^van der Pol equation DAE equivalent
  8:    Processors: 1
  9: */
 10: /* ------------------------------------------------------------------------

 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.

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

 24:   ------------------------------------------------------------------------- */

 26: #include <petscts.h>

 28: typedef struct _n_User *User;
 29: struct _n_User {
 30:   PetscReal mu;
 31:   PetscReal next_output;
 32: };

 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;

 41:   VecGetArrayRead(X,&x);
 42:   VecGetArrayRead(Xdot,&xdot);
 43:   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:   VecRestoreArrayRead(X,&x);
 47:   VecRestoreArrayRead(Xdot,&xdot);
 48:   VecRestoreArray(F,&f);
 49:   return 0;
 50: }

 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];

 60:   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:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 64:   VecRestoreArrayRead(X,&x);

 66:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 67:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 68:   if (A != B) {
 69:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 70:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 71:   }
 72:   return 0;
 73: }

 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;

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:      Initialize program
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:   PetscInitialize(&argc,&argv,NULL,help);
 91:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:     Set runtime options
 96:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   user.next_output = 0.0;
 98:   user.mu          = 1.0e6;
 99:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
100:   PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);
101:   PetscOptionsEnd();

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:     Create necessary matrix and vectors, solve same ODE on every process
105:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106:   MatCreate(PETSC_COMM_WORLD,&A);
107:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
108:   MatSetFromOptions(A);
109:   MatSetUp(A);

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

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:      Create timestepping solver context
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   TSCreate(PETSC_COMM_WORLD,&ts);
117:   TSSetType(ts,TSBEULER);
118:   TSSetIFunction(ts,NULL,IFunction,&user);
119:   TSSetIJacobian(ts,A,A,IJacobian,&user);

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

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Set runtime options
134:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   TSSetFromOptions(ts);

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Solve nonlinear system
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   TSSolve(ts,x);
141:   TSGetSolveTime(ts,&ftime);
142:   TSGetStepNumber(ts,&steps);
143:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
144:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:      Free work space.  All PETSc objects should be destroyed when they
148:      are no longer needed.
149:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150:   MatDestroy(&A);
151:   VecDestroy(&x);
152:   TSDestroy(&ts);

154:   PetscFinalize();
155:   return(ierr);
156: }

158: /*TEST

160:     build:
161:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) radau5

163:     test:
164:       args: -ts_monitor_solution -ts_type radau5

166: TEST*/