Actual source code: ex20.c

petsc-dev 2014-04-23
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: const PetscReal mu = 1.0e6;

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

 86: /*
 87: *  User-defined routines
 88: */
 91: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 92: {
 94:   User           user = (User)ctx;
 95:   PetscScalar    *x,*f;

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

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

116:   VecGetArray(X,&x);
117:   VecGetArray(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:   VecRestoreArray(X,&x);
122:   VecRestoreArray(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: {
132:   User           user     = (User)ctx;
133:   PetscInt       rowcol[] = {0,1};
134:   PetscScalar    *x,J[2][2];

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

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

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

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

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

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

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

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

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

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

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

228:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229:     Set runtime options
230:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231:   user.imex        = PETSC_TRUE;
232:   user.next_output = 0.0;
233:   user.mu          = 1.0e6;
234:   PetscOptionsGetBool(NULL,"-imex",&user.imex,NULL);
235:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
236:   PetscOptionsReal("-mu","Stiffness parameter","<1.0e6>",user.mu,&user.mu,NULL);

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

246:   MatGetVecs(A,&x,NULL);

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

257:   TSSetDuration(ts,PETSC_DEFAULT,ftime);
258:   if (monitor) {
259:     TSMonitorSet(ts,Monitor,&user,NULL);
260:   }

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

270:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271:      Set runtime options
272:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
273:   TSSetFromOptions(ts);

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

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

292:   PetscFinalize();
293:   return(0);
294: }