Actual source code: ex20.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) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
 17:    and
 18:        \mu = 10^6 ( y'(0) ~ -0.6666665432100101).
 19:    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.

 21:    Notes:
 22:    This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form.

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

 26: #include <petscts.h>

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

 34: /*
 35:    User-defined routines
 36: */
 37: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 38: {
 39:   User              user = (User)ctx;
 40:   PetscScalar       *f;
 41:   const PetscScalar *x;

 44:   VecGetArrayRead(X,&x);
 45:   VecGetArray(F,&f);
 46:   f[0] = x[1];
 47:   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
 48:   VecRestoreArrayRead(X,&x);
 49:   VecRestoreArray(F,&f);
 50:   return 0;
 51: }

 53: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 54: {
 55:   User              user = (User)ctx;
 56:   const PetscScalar *x,*xdot;
 57:   PetscScalar       *f;

 60:   VecGetArrayRead(X,&x);
 61:   VecGetArrayRead(Xdot,&xdot);
 62:   VecGetArray(F,&f);
 63:   f[0] = xdot[0] - x[1];
 64:   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
 65:   VecRestoreArrayRead(X,&x);
 66:   VecRestoreArrayRead(Xdot,&xdot);
 67:   VecRestoreArray(F,&f);
 68:   return 0;
 69: }

 71: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 72: {
 73:   User              user     = (User)ctx;
 74:   PetscInt          rowcol[] = {0,1};
 75:   const PetscScalar *x;
 76:   PetscScalar       J[2][2];

 79:   VecGetArrayRead(X,&x);
 80:   J[0][0] = a;     J[0][1] = -1.0;
 81:   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:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 83:   VecRestoreArrayRead(X,&x);

 85:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 86:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 87:   if (A != B) {
 88:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 89:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 90:   }
 91:   return 0;
 92: }

 94: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
 95: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
 96: {
 97:   PetscErrorCode    ierr;
 98:   const PetscScalar *x;
 99:   PetscReal         tfinal, dt;
100:   User              user = (User)ctx;
101:   Vec               interpolatedX;

104:   TSGetTimeStep(ts,&dt);
105:   TSGetMaxTime(ts,&tfinal);

107:   while (user->next_output <= t && user->next_output <= tfinal) {
108:     VecDuplicate(X,&interpolatedX);
109:     TSInterpolate(ts,user->next_output,interpolatedX);
110:     VecGetArrayRead(interpolatedX,&x);
111:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
112:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
113:                        (double)PetscRealPart(x[1]));
114:     VecRestoreArrayRead(interpolatedX,&x);
115:     VecDestroy(&interpolatedX);
116:     user->next_output += 0.1;
117:   }
118:   return 0;
119: }

121: int main(int argc,char **argv)
122: {
123:   TS             ts;            /* nonlinear solver */
124:   Vec            x;             /* solution, residual vectors */
125:   Mat            A;             /* Jacobian matrix */
126:   PetscInt       steps;
127:   PetscReal      ftime = 0.5;
128:   PetscBool      monitor = PETSC_FALSE,implicitform = PETSC_TRUE;
129:   PetscScalar    *x_ptr;
130:   PetscMPIInt    size;
131:   struct _n_User user;

134:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135:      Initialize program
136:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137:   PetscInitialize(&argc,&argv,NULL,help);
138:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:     Set runtime options
143:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144:   user.next_output = 0.0;
145:   user.mu          = 1.0e3;
146:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
147:   PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);
148:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
149:   PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);
150:   PetscOptionsEnd();

152:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153:     Create necessary matrix and vectors, solve same ODE on every process
154:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155:   MatCreate(PETSC_COMM_WORLD,&A);
156:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
157:   MatSetFromOptions(A);
158:   MatSetUp(A);

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

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Create timestepping solver context
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   TSCreate(PETSC_COMM_WORLD,&ts);
166:   if (implicitform) {
167:     TSSetIFunction(ts,NULL,IFunction,&user);
168:     TSSetIJacobian(ts,A,A,IJacobian,&user);
169:     TSSetType(ts,TSBEULER);
170:   } else {
171:     TSSetRHSFunction(ts,NULL,RHSFunction,&user);
172:     TSSetType(ts,TSRK);
173:   }
174:   TSSetMaxTime(ts,ftime);
175:   TSSetTimeStep(ts,0.001);
176:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
177:   if (monitor) {
178:     TSMonitorSet(ts,Monitor,&user,NULL);
179:   }

181:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182:      Set initial conditions
183:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
184:   VecGetArray(x,&x_ptr);
185:   x_ptr[0] = 2.0;
186:   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
187:   VecRestoreArray(x,&x_ptr);

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190:      Set runtime options
191:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192:   TSSetFromOptions(ts);

194:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195:      Solve nonlinear system
196:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197:   TSSolve(ts,x);
198:   TSGetSolveTime(ts,&ftime);
199:   TSGetStepNumber(ts,&steps);
200:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
201:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

203:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:      Free work space.  All PETSc objects should be destroyed when they
205:      are no longer needed.
206:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207:   MatDestroy(&A);
208:   VecDestroy(&x);
209:   TSDestroy(&ts);

211:   PetscFinalize();
212:   return(ierr);
213: }

215: /*TEST

217:     test:
218:       requires: !single
219:       args: -mu 1e6

221:     test:
222:       requires: !single
223:       suffix: 2
224:       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp

226:     test:
227:       requires: !single
228:       suffix: 3
229:       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312

231: TEST*/