Actual source code: ex4.c

petsc-3.4.5 2014-06-29
  2: static char help[] ="Solves a simple time-dependent linear PDE (the heat equation).\n\
  3: Input parameters include:\n\
  4:   -m <points>, where <points> = number of grid points\n\
  5:   -time_dependent_rhs : Treat the problem as having a time-dependent right-hand side\n\
  6:   -debug              : Activate debugging printouts\n\
  7:   -nox                : Deactivate x-window graphics\n\n";

  9: /*
 10:    Concepts: TS^time-dependent linear problems
 11:    Concepts: TS^heat equation
 12:    Concepts: TS^diffusion equation
 13:    Processors: n
 14: */

 16: /* ------------------------------------------------------------------------

 18:    This program solves the one-dimensional heat equation (also called the
 19:    diffusion equation),
 20:        u_t = u_xx,
 21:    on the domain 0 <= x <= 1, with the boundary conditions
 22:        u(t,0) = 0, u(t,1) = 0,
 23:    and the initial condition
 24:        u(0,x) = sin(6*pi*x) + 3*sin(2*pi*x).
 25:    This is a linear, second-order, parabolic equation.

 27:    We discretize the right-hand side using finite differences with
 28:    uniform grid spacing h:
 29:        u_xx = (u_{i+1} - 2u_{i} + u_{i-1})/(h^2)
 30:    We then demonstrate time evolution using the various TS methods by
 31:    running the program via
 32:        mpiexec -n <procs> ex3 -ts_type <timestepping solver>

 34:    We compare the approximate solution with the exact solution, given by
 35:        u_exact(x,t) = exp(-36*pi*pi*t) * sin(6*pi*x) +
 36:                       3*exp(-4*pi*pi*t) * sin(2*pi*x)

 38:    Notes:
 39:    This code demonstrates the TS solver interface to two variants of
 40:    linear problems, u_t = f(u,t), namely
 41:      - time-dependent f:   f(u,t) is a function of t
 42:      - time-independent f: f(u,t) is simply f(u)

 44:     The uniprocessor version of this code is ts/examples/tutorials/ex3.c

 46:   ------------------------------------------------------------------------- */

 48: /*
 49:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs) to manage
 50:    the parallel grid.  Include "petscts.h" so that we can use TS solvers.
 51:    Note that this file automatically includes:
 52:      petscsys.h       - base PETSc routines   petscvec.h  - vectors
 53:      petscmat.h  - matrices
 54:      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
 55:      petscviewer.h - viewers               petscpc.h   - preconditioners
 56:      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
 57: */

 59: #include <petscdmda.h>
 60: #include <petscts.h>
 61: #include <petscdraw.h>

 63: /*
 64:    User-defined application context - contains data needed by the
 65:    application-provided call-back routines.
 66: */
 67: typedef struct {
 68:   MPI_Comm    comm;              /* communicator */
 69:   DM          da;                /* distributed array data structure */
 70:   Vec         localwork;         /* local ghosted work vector */
 71:   Vec         u_local;           /* local ghosted approximate solution vector */
 72:   Vec         solution;          /* global exact solution vector */
 73:   PetscInt    m;                 /* total number of grid points */
 74:   PetscReal   h;                 /* mesh width h = 1/(m-1) */
 75:   PetscBool   debug;             /* flag (1 indicates activation of debugging printouts) */
 76:   PetscViewer viewer1,viewer2;  /* viewers for the solution and error */
 77:   PetscReal   norm_2,norm_max;  /* error norms */
 78: } AppCtx;

 80: /*
 81:    User-defined routines
 82: */
 83: extern PetscErrorCode InitialConditions(Vec,AppCtx*);
 84: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
 85: extern PetscErrorCode RHSFunctionHeat(TS,PetscReal,Vec,Vec,void*);
 86: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
 87: extern PetscErrorCode ExactSolution(PetscReal,Vec,AppCtx*);

 91: int main(int argc,char **argv)
 92: {
 93:   AppCtx         appctx;                 /* user-defined application context */
 94:   TS             ts;                     /* timestepping context */
 95:   Mat            A;                      /* matrix data structure */
 96:   Vec            u;                      /* approximate solution vector */
 97:   PetscReal      time_total_max = 1.0;   /* default max total time */
 98:   PetscInt       time_steps_max = 100;   /* default max timesteps */
 99:   PetscDraw      draw;                   /* drawing context */
101:   PetscInt       steps,m;
102:   PetscMPIInt    size;
103:   PetscReal      dt,ftime;
104:   PetscBool      flg;
105:   TSProblemType  tsproblem = TS_LINEAR;

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:      Initialize program and set problem parameters
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

111:   PetscInitialize(&argc,&argv,(char*)0,help);
112:   appctx.comm = PETSC_COMM_WORLD;

114:   m               = 60;
115:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
116:   PetscOptionsHasName(NULL,"-debug",&appctx.debug);
117:   appctx.m        = m;
118:   appctx.h        = 1.0/(m-1.0);
119:   appctx.norm_2   = 0.0;
120:   appctx.norm_max = 0.0;

122:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
123:   PetscPrintf(PETSC_COMM_WORLD,"Solving a linear TS problem, number of processors = %d\n",size);

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126:      Create vector data structures
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128:   /*
129:      Create distributed array (DMDA) to manage parallel grid and vectors
130:      and to set up the ghost point communication pattern.  There are M
131:      total grid values spread equally among all the processors.
132:   */

134:   DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,m,1,1,NULL,&appctx.da);

136:   /*
137:      Extract global and local vectors from DMDA; we use these to store the
138:      approximate solution.  Then duplicate these for remaining vectors that
139:      have the same types.
140:   */
141:   DMCreateGlobalVector(appctx.da,&u);
142:   DMCreateLocalVector(appctx.da,&appctx.u_local);

144:   /*
145:      Create local work vector for use in evaluating right-hand-side function;
146:      create global work vector for storing exact solution.
147:   */
148:   VecDuplicate(appctx.u_local,&appctx.localwork);
149:   VecDuplicate(u,&appctx.solution);

151:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152:      Set up displays to show graphs of the solution and error
153:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

155:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,380,400,160,&appctx.viewer1);
156:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
157:   PetscDrawSetDoubleBuffer(draw);
158:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,0,400,160,&appctx.viewer2);
159:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
160:   PetscDrawSetDoubleBuffer(draw);

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Create timestepping solver context
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

166:   TSCreate(PETSC_COMM_WORLD,&ts);

168:   flg  = PETSC_FALSE;
169:   PetscOptionsGetBool(NULL,"-nonlinear",&flg,NULL);
170:   TSSetProblemType(ts,flg ? TS_NONLINEAR : TS_LINEAR);

172:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173:      Set optional user-defined monitoring routine
174:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175:   TSMonitorSet(ts,Monitor,&appctx,NULL);

177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

179:      Create matrix data structure; set matrix evaluation routine.
180:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

182:   MatCreate(PETSC_COMM_WORLD,&A);
183:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);
184:   MatSetFromOptions(A);
185:   MatSetUp(A);

187:   flg  = PETSC_FALSE;
188:   PetscOptionsGetBool(NULL,"-time_dependent_rhs",&flg,NULL);
189:   if (flg) {
190:     /*
191:        For linear problems with a time-dependent f(u,t) in the equation
192:        u_t = f(u,t), the user provides the discretized right-hand-side
193:        as a time-dependent matrix.
194:     */
195:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
196:     TSSetRHSJacobian(ts,A,A,RHSMatrixHeat,&appctx);
197:   } else {
198:     /*
199:        For linear problems with a time-independent f(u) in the equation
200:        u_t = f(u), the user provides the discretized right-hand-side
201:        as a matrix only once, and then sets a null matrix evaluation
202:        routine.
203:     */
204:     MatStructure A_structure;
205:     RHSMatrixHeat(ts,0.0,u,&A,&A,&A_structure,&appctx);
206:     TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
207:     TSSetRHSJacobian(ts,A,A,TSComputeRHSJacobianConstant,&appctx);
208:   }

210:   if (tsproblem == TS_NONLINEAR) {
211:     SNES snes;
212:     TSSetRHSFunction(ts,NULL,RHSFunctionHeat,&appctx);
213:     TSGetSNES(ts,&snes);
214:     SNESSetJacobian(snes,NULL,NULL,SNESComputeJacobianDefault,NULL);
215:   }

217:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218:      Set solution vector and initial timestep
219:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

221:   dt   = appctx.h*appctx.h/2.0;
222:   TSSetInitialTimeStep(ts,0.0,dt);
223:   TSSetSolution(ts,u);

225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      Customize timestepping solver:
227:        - Set the solution method to be the Backward Euler method.
228:        - Set timestepping duration info
229:      Then set runtime options, which can override these defaults.
230:      For example,
231:           -ts_max_steps <maxsteps> -ts_final_time <maxtime>
232:      to override the defaults set by TSSetDuration().
233:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

235:   TSSetDuration(ts,time_steps_max,time_total_max);
236:   TSSetFromOptions(ts);

238:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239:      Solve the problem
240:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

242:   /*
243:      Evaluate initial conditions
244:   */
245:   InitialConditions(u,&appctx);

247:   /*
248:      Run the timestepping solver
249:   */
250:   TSSolve(ts,u);
251:   TSGetSolveTime(ts,&ftime);
252:   TSGetTimeStepNumber(ts,&steps);

254:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255:      View timestepping solver info
256:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
257:   PetscPrintf(PETSC_COMM_WORLD,"Total timesteps %D, Final time %G\n",steps,ftime);
258:   PetscPrintf(PETSC_COMM_WORLD,"Avg. error (2 norm) = %G Avg. error (max norm) = %G\n",appctx.norm_2/steps,appctx.norm_max/steps);

260:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261:      Free work space.  All PETSc objects should be destroyed when they
262:      are no longer needed.
263:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

265:   TSDestroy(&ts);
266:   MatDestroy(&A);
267:   VecDestroy(&u);
268:   PetscViewerDestroy(&appctx.viewer1);
269:   PetscViewerDestroy(&appctx.viewer2);
270:   VecDestroy(&appctx.localwork);
271:   VecDestroy(&appctx.solution);
272:   VecDestroy(&appctx.u_local);
273:   DMDestroy(&appctx.da);

275:   /*
276:      Always call PetscFinalize() before exiting a program.  This routine
277:        - finalizes the PETSc libraries as well as MPI
278:        - provides summary and diagnostic information if certain runtime
279:          options are chosen (e.g., -log_summary).
280:   */
281:   PetscFinalize();
282:   return 0;
283: }
284: /* --------------------------------------------------------------------- */
287: /*
288:    InitialConditions - Computes the solution at the initial time.

290:    Input Parameter:
291:    u - uninitialized solution vector (global)
292:    appctx - user-defined application context

294:    Output Parameter:
295:    u - vector with solution at initial time (global)
296: */
297: PetscErrorCode InitialConditions(Vec u,AppCtx *appctx)
298: {
299:   PetscScalar    *u_localptr,h = appctx->h;
300:   PetscInt       i,mybase,myend;

303:   /*
304:      Determine starting point of each processor's range of
305:      grid values.
306:   */
307:   VecGetOwnershipRange(u,&mybase,&myend);

309:   /*
310:     Get a pointer to vector data.
311:     - For default PETSc vectors, VecGetArray() returns a pointer to
312:       the data array.  Otherwise, the routine is implementation dependent.
313:     - You MUST call VecRestoreArray() when you no longer need access to
314:       the array.
315:     - Note that the Fortran interface to VecGetArray() differs from the
316:       C version.  See the users manual for details.
317:   */
318:   VecGetArray(u,&u_localptr);

320:   /*
321:      We initialize the solution array by simply writing the solution
322:      directly into the array locations.  Alternatively, we could use
323:      VecSetValues() or VecSetValuesLocal().
324:   */
325:   for (i=mybase; i<myend; i++) u_localptr[i-mybase] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);

327:   /*
328:      Restore vector
329:   */
330:   VecRestoreArray(u,&u_localptr);

332:   /*
333:      Print debugging information if desired
334:   */
335:   if (appctx->debug) {
336:     PetscPrintf(appctx->comm,"initial guess vector\n");
337:     VecView(u,PETSC_VIEWER_STDOUT_WORLD);
338:   }

340:   return 0;
341: }
342: /* --------------------------------------------------------------------- */
345: /*
346:    ExactSolution - Computes the exact solution at a given time.

348:    Input Parameters:
349:    t - current time
350:    solution - vector in which exact solution will be computed
351:    appctx - user-defined application context

353:    Output Parameter:
354:    solution - vector with the newly computed exact solution
355: */
356: PetscErrorCode ExactSolution(PetscReal t,Vec solution,AppCtx *appctx)
357: {
358:   PetscScalar    *s_localptr,h = appctx->h,ex1,ex2,sc1,sc2;
359:   PetscInt       i,mybase,myend;

362:   /*
363:      Determine starting and ending points of each processor's
364:      range of grid values
365:   */
366:   VecGetOwnershipRange(solution,&mybase,&myend);

368:   /*
369:      Get a pointer to vector data.
370:   */
371:   VecGetArray(solution,&s_localptr);

373:   /*
374:      Simply write the solution directly into the array locations.
375:      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
376:   */
377:   ex1 = exp(-36.*PETSC_PI*PETSC_PI*t); ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
378:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
379:   for (i=mybase; i<myend; i++) s_localptr[i-mybase] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;

381:   /*
382:      Restore vector
383:   */
384:   VecRestoreArray(solution,&s_localptr);
385:   return 0;
386: }
387: /* --------------------------------------------------------------------- */
390: /*
391:    Monitor - User-provided routine to monitor the solution computed at
392:    each timestep.  This example plots the solution and computes the
393:    error in two different norms.

395:    Input Parameters:
396:    ts     - the timestep context
397:    step   - the count of the current step (with 0 meaning the
398:              initial condition)
399:    time   - the current time
400:    u      - the solution at this timestep
401:    ctx    - the user-provided context for this monitoring routine.
402:             In this case we use the application context which contains
403:             information about the problem size, workspace and the exact
404:             solution.
405: */
406: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx)
407: {
408:   AppCtx         *appctx = (AppCtx*) ctx;   /* user-defined application context */
410:   PetscReal      norm_2,norm_max;

412:   /*
413:      View a graph of the current iterate
414:   */
415:   VecView(u,appctx->viewer2);

417:   /*
418:      Compute the exact solution
419:   */
420:   ExactSolution(time,appctx->solution,appctx);

422:   /*
423:      Print debugging information if desired
424:   */
425:   if (appctx->debug) {
426:     PetscPrintf(appctx->comm,"Computed solution vector\n");
427:     VecView(u,PETSC_VIEWER_STDOUT_WORLD);
428:     PetscPrintf(appctx->comm,"Exact solution vector\n");
429:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);
430:   }

432:   /*
433:      Compute the 2-norm and max-norm of the error
434:   */
435:   VecAXPY(appctx->solution,-1.0,u);
436:   VecNorm(appctx->solution,NORM_2,&norm_2);
437:   norm_2 = PetscSqrtReal(appctx->h)*norm_2;
438:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

440:   /*
441:      PetscPrintf() causes only the first processor in this
442:      communicator to print the timestep information.
443:   */
444:   PetscPrintf(appctx->comm,"Timestep %D: time = %G 2-norm error = %G max norm error = %G\n",
445:                      step,time,norm_2,norm_max);
446:   appctx->norm_2   += norm_2;
447:   appctx->norm_max += norm_max;

449:   /*
450:      View a graph of the error
451:   */
452:   VecView(appctx->solution,appctx->viewer1);

454:   /*
455:      Print debugging information if desired
456:   */
457:   if (appctx->debug) {
458:     PetscPrintf(appctx->comm,"Error vector\n");
459:     VecView(appctx->solution,PETSC_VIEWER_STDOUT_WORLD);
460:   }

462:   return 0;
463: }

465: /* --------------------------------------------------------------------- */
468: /*
469:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
470:    matrix for the heat equation.

472:    Input Parameters:
473:    ts - the TS context
474:    t - current time
475:    global_in - global input vector
476:    dummy - optional user-defined context, as set by TSetRHSJacobian()

478:    Output Parameters:
479:    AA - Jacobian matrix
480:    BB - optionally different preconditioning matrix
481:    str - flag indicating matrix structure

483:   Notes:
484:   RHSMatrixHeat computes entries for the locally owned part of the system.
485:    - Currently, all PETSc parallel matrix formats are partitioned by
486:      contiguous chunks of rows across the processors.
487:    - Each processor needs to insert only elements that it owns
488:      locally (but any non-local elements will be sent to the
489:      appropriate processor during matrix assembly).
490:    - Always specify global row and columns of matrix entries when
491:      using MatSetValues(); we could alternatively use MatSetValuesLocal().
492:    - Here, we set all entries for a particular row at once.
493:    - Note that MatSetValues() uses 0-based row and column numbers
494:      in Fortran as well as in C.
495: */
496: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec X,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
497: {
498:   Mat            A       = *AA;              /* Jacobian matrix */
499:   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
501:   PetscInt       i,mstart,mend,idx[3];
502:   PetscScalar    v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

504:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
505:      Compute entries for the locally owned part of the matrix
506:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

508:   MatGetOwnershipRange(A,&mstart,&mend);

510:   /*
511:      Set matrix rows corresponding to boundary data
512:   */

514:   if (mstart == 0) {  /* first processor only */
515:     v[0] = 1.0;
516:     MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
517:     mstart++;
518:   }

520:   if (mend == appctx->m) { /* last processor only */
521:     mend--;
522:     v[0] = 1.0;
523:     MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
524:   }

526:   /*
527:      Set matrix rows corresponding to interior data.  We construct the
528:      matrix one row at a time.
529:   */
530:   v[0] = sone; v[1] = stwo; v[2] = sone;
531:   for (i=mstart; i<mend; i++) {
532:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
533:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
534:   }

536:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
537:      Complete the matrix assembly process and set some options
538:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
539:   /*
540:      Assemble matrix, using the 2-step process:
541:        MatAssemblyBegin(), MatAssemblyEnd()
542:      Computations can be done while messages are in transition
543:      by placing code between these two statements.
544:   */
545:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
546:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

548:   /*
549:      Set flag to indicate that the Jacobian matrix retains an identical
550:      nonzero structure throughout all timestepping iterations (although the
551:      values of the entries change). Thus, we can save some work in setting
552:      up the preconditioner (e.g., no need to redo symbolic factorization for
553:      ILU/ICC preconditioners).
554:       - If the nonzero structure of the matrix is different during
555:         successive linear solves, then the flag DIFFERENT_NONZERO_PATTERN
556:         must be used instead.  If you are unsure whether the matrix
557:         structure has changed or not, use the flag DIFFERENT_NONZERO_PATTERN.
558:       - Caution:  If you specify SAME_NONZERO_PATTERN, PETSc
559:         believes your assertion and does not check the structure
560:         of the matrix.  If you erroneously claim that the structure
561:         is the same when it actually is not, the new preconditioner
562:         will not function correctly.  Thus, use this optimization
563:         feature with caution!
564:   */
565:   *str = SAME_NONZERO_PATTERN;

567:   /*
568:      Set and option to indicate that we will never add a new nonzero location
569:      to the matrix. If we do, it will generate an error.
570:   */
571:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

573:   return 0;
574: }

578: PetscErrorCode RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
579: {
581:   Mat            A;
582:   MatStructure   A_structure;

585:   TSGetRHSJacobian(ts,&A,NULL,NULL,&ctx);
586:   RHSMatrixHeat(ts,t,globalin,&A,NULL,&A_structure,ctx);
587:   /* MatView(A,PETSC_VIEWER_STDOUT_WORLD); */
588:   MatMult(A,globalin,globalout);
589:   return(0);
590: }