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:   PetscErrorCode    ierr;
 40:   User              user = (User)ctx;
 41:   PetscScalar       *f;
 42:   const PetscScalar *x;

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

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

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

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

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

 88:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 89:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 90:   if (A != B) {
 91:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 92:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 93:   }
 94:   return(0);
 95: }

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

107:   TSGetTimeStep(ts,&dt);
108:   TSGetMaxTime(ts,&tfinal);

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

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

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Initialize program
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
141:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
142:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");

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

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

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

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

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

192:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193:      Set runtime options
194:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195:   TSSetFromOptions(ts);

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

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

214:   PetscFinalize();
215:   return(ierr);
216: }

218: /*TEST

220:     test:
221:       requires: !single
222:       args: -mu 1e6

224:     test:
225:       requires: !single
226:       suffix: 2
227:       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp

229:     test:
230:       requires: !single
231:       suffix: 3
232:       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312

234: TEST*/