Actual source code: ex20.c

petsc-3.7.3 2016-08-01
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) = -6.666665432100101e-01,
 17:    and
 18:        mu = 10^6.
 19:    This is a nonlinear equation.

 21:    Notes:
 22:    This code demonstrates the TS solver interface to a variant of
 23:    linear problems, u_t = f(u,t), namely turning (1) into a system of
 24:    first order differential equations,

 26:    [ y' ] = [          z          ]
 27:    [ z' ]   [     mu[(1-y^2)z-y]  ]

 29:    which then we can write as a vector equation

 31:    [ u_1' ] = [      u_2              ]  (2)
 32:    [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1]  ]

 34:    which is now in the desired form of u_t = f(u,t). One way that we
 35:    can split f(u,t) in (2) is to split by component,

 37:    [ u_1' ] = [  u_2 ] + [       0              ]
 38:    [ u_2' ]   [  0   ]   [ mu[(1-u_1^2)u_2-u_1] ]

 40:    where

 42:    [ F(u,t) ] = [  u_2 ]
 43:                 [  0   ]

 45:    and

 47:    [ G(u',u,t) ] = [ u_1' ] - [            0         ]
 48:                    [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1] ]

 50:    Using the definition of the Jacobian of G (from the PETSc user manual),
 51:    in the equation G(u',u,t) = F(u,t),

 53:               dG   dG
 54:    J(G) = a * -- - --
 55:               du'  du

 57:    where d is the partial derivative. In this example,

 59:    dG   [ 1 ; 0 ]
 60:    -- = [       ]
 61:    du'  [ 0 ; 1 ]

 63:    dG   [ 0                       ;         0         ]
 64:    -- = [                                             ]
 65:    du   [ -mu*(1.0 + 2.0*u_1*u_2) ; mu*(1-u_1*u_1)    ]

 67:    Hence,

 69:           [      a                 ;         0          ]
 70:    J(G) = [                                             ]
 71:           [ mu*(1.0 + 2.0*u_1*u_2) ; a - mu*(1-u_1*u_1) ]

 73:   ------------------------------------------------------------------------- */

 75: #include <petscts.h>

 77: typedef struct _n_User *User;
 78: struct _n_User {
 79:   PetscReal mu;
 80:   PetscBool imex;
 81:   PetscReal next_output;
 82: };

 84: /*
 85: *  User-defined routines
 86: */
 89: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 90: {
 91:   PetscErrorCode    ierr;
 92:   User              user = (User)ctx;
 93:   PetscScalar       *f;
 94:   const PetscScalar *x;

 97:   VecGetArrayRead(X,&x);
 98:   VecGetArray(F,&f);
 99:   f[0] = (user->imex ? x[1] : 0.0);
100:   f[1] = 0.0;
101:   VecRestoreArrayRead(X,&x);
102:   VecRestoreArray(F,&f);
103:   return(0);
104: }

108: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
109: {
110:   PetscErrorCode    ierr;
111:   User              user = (User)ctx;
112:   const PetscScalar *x,*xdot;
113:   PetscScalar       *f;

116:   VecGetArrayRead(X,&x);
117:   VecGetArrayRead(Xdot,&xdot);
118:   VecGetArray(F,&f);
119:   f[0] = xdot[0] - (user->imex ? 0 : x[1]);
120:   f[1] = xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]);
121:   VecRestoreArrayRead(X,&x);
122:   VecRestoreArrayRead(Xdot,&xdot);
123:   VecRestoreArray(F,&f);
124:   return(0);
125: }

129: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
130: {
131:   PetscErrorCode    ierr;
132:   User              user     = (User)ctx;
133:   PetscInt          rowcol[] = {0,1};
134:   const PetscScalar *x;
135:   PetscScalar       J[2][2];

138:   VecGetArrayRead(X,&x);
139:   J[0][0] = a;     J[0][1] = (user->imex ? 0 : -1.0);
140:   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]);
141:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
142:   VecRestoreArrayRead(X,&x);

144:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
145:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
146:   if (A != B) {
147:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
148:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
149:   }
150:   return(0);
151: }

155: /* This is an example of registering an user-provided ARKIMEX scheme */
156: static PetscErrorCode RegisterMyARK2(void)
157: {

161:   {
162:     const PetscReal
163:       A[3][3] = {{0,0,0},
164:                  {0.41421356237309504880,0,0},
165:                  {0.75,0.25,0}},
166:       At[3][3] = {{0,0,0},
167:                   {0.12132034355964257320,0.29289321881345247560,0},
168:                   {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}};
169:     TSARKIMEXRegister("myark2",2,3,&At[0][0],NULL,NULL,&A[0][0],NULL,NULL,NULL,NULL,0,NULL,NULL);
170:   }
171:   return(0);
172: }

176: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
177: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
178: {
179:   PetscErrorCode    ierr;
180:   const PetscScalar *x;
181:   PetscReal         tfinal, dt;
182:   User              user = (User)ctx;
183:   Vec               interpolatedX;

186:   TSGetTimeStep(ts,&dt);
187:   TSGetDuration(ts,NULL,&tfinal);

189:   while (user->next_output <= t && user->next_output <= tfinal) {
190:     VecDuplicate(X,&interpolatedX);
191:     TSInterpolate(ts,user->next_output,interpolatedX);
192:     VecGetArrayRead(interpolatedX,&x);
193:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
194:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
195:                        (double)PetscRealPart(x[1]));
196:     VecRestoreArrayRead(interpolatedX,&x);
197:     VecDestroy(&interpolatedX);
198:     user->next_output += 0.1;
199:   }
200:   return(0);
201: }

205: int main(int argc,char **argv)
206: {
207:   TS             ts;            /* nonlinear solver */
208:   Vec            x;             /* solution, residual vectors */
209:   Mat            A;             /* Jacobian matrix */
210:   PetscInt       steps;
211:   PetscReal      ftime   = 0.5;
212:   PetscBool      monitor = PETSC_FALSE;
213:   PetscScalar    *x_ptr;
214:   PetscMPIInt    size;
215:   struct _n_User user;

218:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219:      Initialize program
220:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
221:   PetscInitialize(&argc,&argv,NULL,help);

223:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
224:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

226:   /* Register user-specified ARKIMEX method */
227:   RegisterMyARK2();

229:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230:     Set runtime options
231:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232:   user.imex        = PETSC_TRUE;
233:   user.next_output = 0.0;
234:   user.mu          = 1.0e6;
235:   PetscOptionsGetBool(NULL,NULL,"-imex",&user.imex,NULL);
236:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
237:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
238:   PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,PETSC_NULL);
239:   PetscOptionsEnd();

241:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242:     Create necessary matrix and vectors, solve same ODE on every process
243:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
244:   MatCreate(PETSC_COMM_WORLD,&A);
245:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
246:   MatSetFromOptions(A);
247:   MatSetUp(A);

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

251:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252:      Create timestepping solver context
253:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
254:   TSCreate(PETSC_COMM_WORLD,&ts);
255:   TSSetType(ts,TSBEULER);
256:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
257:   TSSetIFunction(ts,NULL,IFunction,&user);
258:   TSSetIJacobian(ts,A,A,IJacobian,&user);

260:   TSSetDuration(ts,PETSC_DEFAULT,ftime);
261:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
262:   if (monitor) {
263:     TSMonitorSet(ts,Monitor,&user,NULL);
264:   }

266:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267:      Set initial conditions
268:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269:   VecGetArray(x,&x_ptr);
270:   x_ptr[0] = 2.0;   x_ptr[1] = -6.666665432100101e-01;
271:   VecRestoreArray(x,&x_ptr);
272:   TSSetInitialTimeStep(ts,0.0,.001);

274:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275:      Set runtime options
276:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
277:   TSSetFromOptions(ts);

279:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280:      Solve nonlinear system
281:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
282:   TSSolve(ts,x);
283:   TSGetSolveTime(ts,&ftime);
284:   TSGetTimeStepNumber(ts,&steps);
285:   PetscPrintf(PETSC_COMM_WORLD,"steps %D, ftime %g\n",steps,(double)ftime);
286:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

288:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
289:      Free work space.  All PETSc objects should be destroyed when they
290:      are no longer needed.
291:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
292:   MatDestroy(&A);
293:   VecDestroy(&x);
294:   TSDestroy(&ts);

296:   PetscFinalize();
297:   return(0);
298: }