Actual source code: ex20.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  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: }


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

219: /*TEST

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

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

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

235: TEST*/