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

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

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

 21: #include <petscts.h>

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

 29: /*
 30:    User-defined routines
 31: */
 32: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 33: {
 34:   User              user = (User)ctx;
 35:   PetscScalar       *f;
 36:   const PetscScalar *x;

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

 48: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 49: {
 50:   User              user = (User)ctx;
 51:   const PetscScalar *x,*xdot;
 52:   PetscScalar       *f;

 55:   VecGetArrayRead(X,&x);
 56:   VecGetArrayRead(Xdot,&xdot);
 57:   VecGetArray(F,&f);
 58:   f[0] = xdot[0] - x[1];
 59:   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
 60:   VecRestoreArrayRead(X,&x);
 61:   VecRestoreArrayRead(Xdot,&xdot);
 62:   VecRestoreArray(F,&f);
 63:   return 0;
 64: }

 66: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 67: {
 68:   User              user     = (User)ctx;
 69:   PetscInt          rowcol[] = {0,1};
 70:   const PetscScalar *x;
 71:   PetscScalar       J[2][2];

 74:   VecGetArrayRead(X,&x);
 75:   J[0][0] = a;     J[0][1] = -1.0;
 76:   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]);
 77:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 78:   VecRestoreArrayRead(X,&x);

 80:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 81:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 82:   if (A != B) {
 83:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 84:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 85:   }
 86:   return 0;
 87: }

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

 99:   TSGetTimeStep(ts,&dt);
100:   TSGetMaxTime(ts,&tfinal);

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

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

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130:      Initialize program
131:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132:   PetscInitialize(&argc,&argv,NULL,help);
133:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:     Create necessary matrix and vectors, solve same ODE on every process
149:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150:   MatCreate(PETSC_COMM_WORLD,&A);
151:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
152:   MatSetFromOptions(A);
153:   MatSetUp(A);

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

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

176:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177:      Set initial conditions
178:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
179:   VecGetArray(x,&x_ptr);
180:   x_ptr[0] = 2.0;
181:   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
182:   VecRestoreArray(x,&x_ptr);

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:      Set runtime options
186:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187:   TSSetFromOptions(ts);

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190:      Solve nonlinear system
191:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192:   TSSolve(ts,x);
193:   TSGetSolveTime(ts,&ftime);
194:   TSGetStepNumber(ts,&steps);
195:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
196:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:      Free work space.  All PETSc objects should be destroyed when they
200:      are no longer needed.
201:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202:   MatDestroy(&A);
203:   VecDestroy(&x);
204:   TSDestroy(&ts);

206:   PetscFinalize();
207:   return(ierr);
208: }

210: /*TEST

212:     test:
213:       requires: !single
214:       args: -mu 1e6

216:     test:
217:       requires: !single
218:       suffix: 2
219:       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp

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

226: TEST*/