Actual source code: chwirut2.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: /*
  2:    Include "petsctao.h" so that we can use TAO solvers.  Note that this
  3:    file automatically includes libraries such as:
  4:      petsc.h       - base PETSc routines   petscvec.h - vectors
  5:      petscsys.h    - sysem routines        petscmat.h - matrices
  6:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
  7:      petscviewer.h - viewers               petscpc.h  - preconditioners

  9: */

 11: #include <petsctao.h>
 12: #include <mpi.h>


 15: /*
 16: Description:   These data are the result of a NIST study involving
 17:                ultrasonic calibration.  The response variable is
 18:                ultrasonic response, and the predictor variable is
 19:                metal distance.

 21: Reference:     Chwirut, D., NIST (197?).
 22:                Ultrasonic Reference Block Study.
 23: */

 25: static char help[]="Finds the nonlinear least-squares solution to the model \n\
 26:             y = exp[-b1*x]/(b2+b3*x)  +  e \n";

 28: /* T
 29:    Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
 30:    Routines: TaoCreate();
 31:    Routines: TaoSetType();
 32:    Routines: TaoSetSeparableObjectiveRoutine();
 33:    Routines: TaoSetMonitor();
 34:    Routines: TaoSetInitialVector();
 35:    Routines: TaoSetFromOptions();
 36:    Routines: TaoSolve();
 37:    Routines: TaoDestroy();
 38:    Processors: n
 39: T*/

 41: #define NOBSERVATIONS 214
 42: #define NPARAMETERS 3

 44: #define DIE_TAG 2000
 45: #define IDLE_TAG 1000

 47: /* User-defined application context */
 48: typedef struct {
 49:   /* Working space */
 50:   PetscReal   t[NOBSERVATIONS];   /* array of independent variables of observation */
 51:   PetscReal   y[NOBSERVATIONS];   /* array of dependent variables */
 52:   PetscMPIInt size,rank;
 53: } AppCtx;

 55: /* User provided Routines */
 56: PetscErrorCode InitializeData(AppCtx *user);
 57: PetscErrorCode FormStartingPoint(Vec);
 58: PetscErrorCode EvaluateFunction(Tao, Vec, Vec, void *);
 59: PetscErrorCode MyMonitor(Tao, void*);
 60: PetscErrorCode TaskWorker(AppCtx *user);
 61: PetscErrorCode StopWorkers(AppCtx *user);
 62: PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal*f, AppCtx *user);

 64: /*--------------------------------------------------------------------*/
 67: int main(int argc,char **argv)
 68: {
 70:   Vec            x, f;               /* solution, function */
 71:   Tao            tao;                /* Tao solver context */
 72:   AppCtx         user;               /* user-defined work context */

 74:    /* Initialize TAO and PETSc */
 75:   PetscInitialize(&argc,&argv,(char *)0,help);

 77:   MPI_Comm_size(MPI_COMM_WORLD,&user.size);
 78:   MPI_Comm_rank(MPI_COMM_WORLD,&user.rank);
 79:   InitializeData(&user);

 81:   /* Run optimization on rank 0 */
 82:   if (user.rank == 0) {
 83:     /* Allocate vectors */
 84:     VecCreateSeq(PETSC_COMM_SELF,NPARAMETERS,&x);
 85:     VecCreateSeq(PETSC_COMM_SELF,NOBSERVATIONS,&f);

 87:     /* TAO code begins here */

 89:     /* Create TAO solver and set desired solution method */
 90:     TaoCreate(PETSC_COMM_SELF,&tao);
 91:     TaoSetType(tao,TAOPOUNDERS);

 93:     /* Set the function and Jacobian routines. */
 94:     FormStartingPoint(x);
 95:     TaoSetInitialVector(tao,x);
 96:     TaoSetSeparableObjectiveRoutine(tao,f,EvaluateFunction,(void*)&user);
 97:     TaoSetMonitor(tao,MyMonitor,&user,NULL);


100:     /* Check for any TAO command line arguments */
101:     TaoSetFromOptions(tao);

103:     /* Perform the Solve */
104:     TaoSolve(tao);

106:     /* Free TAO data structures */
107:     TaoDestroy(&tao);

109:     /* Free PETSc data structures */
110:     VecDestroy(&x);
111:     VecDestroy(&f);
112:     StopWorkers(&user);
113:   } else {
114:     TaskWorker(&user);
115:   }
116:   PetscFinalize();
117:   return 0;
118: }

120: /*--------------------------------------------------------------------*/
123: PetscErrorCode EvaluateFunction(Tao tao, Vec X, Vec F, void *ptr)
124: {
125:   AppCtx         *user = (AppCtx *)ptr;
126:   PetscInt       i;
127:   PetscReal      *x,*f;

131:   VecGetArray(X,&x);
132:   VecGetArray(F,&f);
133:   if (user->size == 1) {
134:     /* Single processor */
135:     for (i=0;i<NOBSERVATIONS;i++) {
136:       RunSimulation(x,i,&f[i],user);
137:     }
138:   } else {
139:     /* Multiprocessor master */
140:     PetscMPIInt tag;
141:     PetscInt    finishedtasks,next_task,checkedin;
142:     PetscReal   f_i;
143:     MPI_Status  status;

145:     next_task=0;
146:     finishedtasks=0;
147:     checkedin=0;

149:     while(finishedtasks < NOBSERVATIONS || checkedin < user->size-1) {
150:       MPI_Recv(&f_i,1,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PETSC_COMM_WORLD,&status);
151:       if (status.MPI_TAG == IDLE_TAG) {
152:         checkedin++;
153:       } else {

155:         tag = status.MPI_TAG;
156:         f[tag] = (PetscReal)f_i;
157:         finishedtasks++;
158:       }

160:       if (next_task<NOBSERVATIONS) {
161:         MPI_Send(x,NPARAMETERS,MPIU_REAL,status.MPI_SOURCE,next_task,PETSC_COMM_WORLD);
162:         next_task++;

164:       } else {
165:         /* Send idle message */
166:         MPI_Send(x,NPARAMETERS,MPIU_REAL,status.MPI_SOURCE,IDLE_TAG,PETSC_COMM_WORLD);
167:       }
168:     }
169:   }
170:   VecRestoreArray(X,&x);
171:   VecRestoreArray(F,&f);
172:   PetscLogFlops(6*NOBSERVATIONS);
173:   return(0);
174: }

176: /* ------------------------------------------------------------ */
179: PetscErrorCode FormStartingPoint(Vec X)
180: {
181:   PetscReal      *x;

185:   VecGetArray(X,&x);
186:   x[0] = 0.15;
187:   x[1] = 0.008;
188:   x[2] = 0.010;
189:   VecRestoreArray(X,&x);
190:   return(0);
191: }

193: /* ---------------------------------------------------------------------- */
196: PetscErrorCode InitializeData(AppCtx *user)
197: {
198:   PetscReal *t=user->t,*y=user->y;
199:   PetscInt  i=0;

202:   y[i] =   92.9000;   t[i++] =  0.5000;
203:   y[i] =    78.7000;  t[i++] =   0.6250;
204:   y[i] =    64.2000;  t[i++] =   0.7500;
205:   y[i] =    64.9000;  t[i++] =   0.8750;
206:   y[i] =    57.1000;  t[i++] =   1.0000;
207:   y[i] =    43.3000;  t[i++] =   1.2500;
208:   y[i] =    31.1000;   t[i++] =  1.7500;
209:   y[i] =    23.6000;   t[i++] =  2.2500;
210:   y[i] =    31.0500;   t[i++] =  1.7500;
211:   y[i] =    23.7750;   t[i++] =  2.2500;
212:   y[i] =    17.7375;   t[i++] =  2.7500;
213:   y[i] =    13.8000;   t[i++] =  3.2500;
214:   y[i] =    11.5875;   t[i++] =  3.7500;
215:   y[i] =     9.4125;   t[i++] =  4.2500;
216:   y[i] =     7.7250;   t[i++] =  4.7500;
217:   y[i] =     7.3500;   t[i++] =  5.2500;
218:   y[i] =     8.0250;   t[i++] =  5.7500;
219:   y[i] =    90.6000;   t[i++] =  0.5000;
220:   y[i] =    76.9000;   t[i++] =  0.6250;
221:   y[i] =    71.6000;   t[i++] = 0.7500;
222:   y[i] =    63.6000;   t[i++] =  0.8750;
223:   y[i] =    54.0000;   t[i++] =  1.0000;
224:   y[i] =    39.2000;   t[i++] =  1.2500;
225:   y[i] =    29.3000;   t[i++] = 1.7500;
226:   y[i] =    21.4000;   t[i++] =  2.2500;
227:   y[i] =    29.1750;   t[i++] =  1.7500;
228:   y[i] =    22.1250;   t[i++] =  2.2500;
229:   y[i] =    17.5125;   t[i++] =  2.7500;
230:   y[i] =    14.2500;   t[i++] =  3.2500;
231:   y[i] =     9.4500;   t[i++] =  3.7500;
232:   y[i] =     9.1500;   t[i++] =  4.2500;
233:   y[i] =     7.9125;   t[i++] =  4.7500;
234:   y[i] =     8.4750;   t[i++] =  5.2500;
235:   y[i] =     6.1125;   t[i++] =  5.7500;
236:   y[i] =    80.0000;   t[i++] =  0.5000;
237:   y[i] =    79.0000;   t[i++] =  0.6250;
238:   y[i] =    63.8000;   t[i++] =  0.7500;
239:   y[i] =    57.2000;   t[i++] =  0.8750;
240:   y[i] =    53.2000;   t[i++] =  1.0000;
241:   y[i] =   42.5000;   t[i++] =  1.2500;
242:   y[i] =   26.8000;   t[i++] =  1.7500;
243:   y[i] =    20.4000;   t[i++] =  2.2500;
244:   y[i] =    26.8500;  t[i++] =   1.7500;
245:   y[i] =    21.0000;  t[i++] =   2.2500;
246:   y[i] =    16.4625;  t[i++] =   2.7500;
247:   y[i] =    12.5250;  t[i++] =   3.2500;
248:   y[i] =    10.5375;  t[i++] =   3.7500;
249:   y[i] =     8.5875;  t[i++] =   4.2500;
250:   y[i] =     7.1250;  t[i++] =   4.7500;
251:   y[i] =     6.1125;  t[i++] =   5.2500;
252:   y[i] =     5.9625;  t[i++] =   5.7500;
253:   y[i] =    74.1000;  t[i++] =   0.5000;
254:   y[i] =    67.3000;  t[i++] =   0.6250;
255:   y[i] =    60.8000;  t[i++] =   0.7500;
256:   y[i] =    55.5000;  t[i++] =   0.8750;
257:   y[i] =    50.3000;  t[i++] =   1.0000;
258:   y[i] =    41.0000;  t[i++] =   1.2500;
259:   y[i] =    29.4000;  t[i++] =   1.7500;
260:   y[i] =    20.4000;  t[i++] =   2.2500;
261:   y[i] =    29.3625;  t[i++] =   1.7500;
262:   y[i] =    21.1500;  t[i++] =   2.2500;
263:   y[i] =    16.7625;  t[i++] =   2.7500;
264:   y[i] =    13.2000;  t[i++] =   3.2500;
265:   y[i] =    10.8750;  t[i++] =   3.7500;
266:   y[i] =     8.1750;  t[i++] =   4.2500;
267:   y[i] =     7.3500;  t[i++] =   4.7500;
268:   y[i] =     5.9625;  t[i++] =  5.2500;
269:   y[i] =     5.6250;  t[i++] =   5.7500;
270:   y[i] =    81.5000;  t[i++] =    .5000;
271:   y[i] =    62.4000;  t[i++] =    .7500;
272:   y[i] =    32.5000;  t[i++] =   1.5000;
273:   y[i] =    12.4100;  t[i++] =   3.0000;
274:   y[i] =    13.1200;  t[i++] =   3.0000;
275:   y[i] =    15.5600;  t[i++] =   3.0000;
276:   y[i] =     5.6300;  t[i++] =   6.0000;
277:   y[i] =    78.0000;   t[i++] =   .5000;
278:   y[i] =    59.9000;  t[i++] =    .7500;
279:   y[i] =    33.2000;  t[i++] =   1.5000;
280:   y[i] =    13.8400;  t[i++] =   3.0000;
281:   y[i] =    12.7500;  t[i++] =   3.0000;
282:   y[i] =    14.6200;  t[i++] =   3.0000;
283:   y[i] =     3.9400;  t[i++] =   6.0000;
284:   y[i] =    76.8000;  t[i++] =    .5000;
285:   y[i] =    61.0000;  t[i++] =    .7500;
286:   y[i] =    32.9000;  t[i++] =   1.5000;
287:   y[i] =   13.8700;   t[i++] = 3.0000;
288:   y[i] =    11.8100;  t[i++] =   3.0000;
289:   y[i] =    13.3100;  t[i++] =   3.0000;
290:   y[i] =     5.4400;  t[i++] =   6.0000;
291:   y[i] =    78.0000;  t[i++] =    .5000;
292:   y[i] =    63.5000;  t[i++] =    .7500;
293:   y[i] =    33.8000;  t[i++] =   1.5000;
294:   y[i] =    12.5600;  t[i++] =   3.0000;
295:   y[i] =     5.6300;  t[i++] =   6.0000;
296:   y[i] =    12.7500;  t[i++] =   3.0000;
297:   y[i] =    13.1200;  t[i++] =   3.0000;
298:   y[i] =     5.4400;  t[i++] =   6.0000;
299:   y[i] =    76.8000;  t[i++] =    .5000;
300:   y[i] =    60.0000;  t[i++] =    .7500;
301:   y[i] =    47.8000;  t[i++] =   1.0000;
302:   y[i] =    32.0000;  t[i++] =   1.5000;
303:   y[i] =    22.2000;  t[i++] =   2.0000;
304:   y[i] =    22.5700;  t[i++] =   2.0000;
305:   y[i] =    18.8200;  t[i++] =   2.5000;
306:   y[i] =    13.9500;  t[i++] =   3.0000;
307:   y[i] =    11.2500;  t[i++] =   4.0000;
308:   y[i] =     9.0000;  t[i++] =   5.0000;
309:   y[i] =     6.6700;  t[i++] =   6.0000;
310:   y[i] =    75.8000;  t[i++] =    .5000;
311:   y[i] =    62.0000;  t[i++] =    .7500;
312:   y[i] =    48.8000;  t[i++] =   1.0000;
313:   y[i] =    35.2000;  t[i++] =   1.5000;
314:   y[i] =    20.0000;  t[i++] =   2.0000;
315:   y[i] =    20.3200;  t[i++] =   2.0000;
316:   y[i] =    19.3100;  t[i++] =   2.5000;
317:   y[i] =    12.7500;  t[i++] =   3.0000;
318:   y[i] =    10.4200;  t[i++] =   4.0000;
319:   y[i] =     7.3100;  t[i++] =   5.0000;
320:   y[i] =     7.4200;  t[i++] =   6.0000;
321:   y[i] =    70.5000;  t[i++] =    .5000;
322:   y[i] =    59.5000;  t[i++] =    .7500;
323:   y[i] =    48.5000;  t[i++] =   1.0000;
324:   y[i] =    35.8000;  t[i++] =   1.5000;
325:   y[i] =    21.0000;  t[i++] =   2.0000;
326:   y[i] =    21.6700;  t[i++] =   2.0000;
327:   y[i] =    21.0000;  t[i++] =   2.5000;
328:   y[i] =    15.6400;  t[i++] =   3.0000;
329:   y[i] =     8.1700;  t[i++] =   4.0000;
330:   y[i] =     8.5500;  t[i++] =   5.0000;
331:   y[i] =    10.1200;  t[i++] =   6.0000;
332:   y[i] =    78.0000;  t[i++] =    .5000;
333:   y[i] =    66.0000;  t[i++] =    .6250;
334:   y[i] =    62.0000;  t[i++] =    .7500;
335:   y[i] =    58.0000;  t[i++] =    .8750;
336:   y[i] =    47.7000;  t[i++] =   1.0000;
337:   y[i] =    37.8000;  t[i++] =   1.2500;
338:   y[i] =    20.2000;  t[i++] =   2.2500;
339:   y[i] =    21.0700;  t[i++] =   2.2500;
340:   y[i] =    13.8700;  t[i++] =   2.7500;
341:   y[i] =     9.6700;  t[i++] =   3.2500;
342:   y[i] =     7.7600;  t[i++] =   3.7500;
343:   y[i] =    5.4400;   t[i++] =  4.2500;
344:   y[i] =    4.8700;   t[i++] =  4.7500;
345:   y[i] =     4.0100;  t[i++] =   5.2500;
346:   y[i] =     3.7500;  t[i++] =   5.7500;
347:   y[i] =    24.1900;  t[i++] =   3.0000;
348:   y[i] =    25.7600;  t[i++] =   3.0000;
349:   y[i] =    18.0700;  t[i++] =   3.0000;
350:   y[i] =    11.8100;  t[i++] =   3.0000;
351:   y[i] =    12.0700;  t[i++] =   3.0000;
352:   y[i] =    16.1200;  t[i++] =   3.0000;
353:   y[i] =    70.8000;  t[i++] =    .5000;
354:   y[i] =    54.7000;  t[i++] =    .7500;
355:   y[i] =    48.0000;  t[i++] =   1.0000;
356:   y[i] =    39.8000;  t[i++] =   1.5000;
357:   y[i] =    29.8000;  t[i++] =   2.0000;
358:   y[i] =    23.7000;  t[i++] =   2.5000;
359:   y[i] =    29.6200;  t[i++] =   2.0000;
360:   y[i] =    23.8100;  t[i++] =   2.5000;
361:   y[i] =    17.7000;  t[i++] =   3.0000;
362:   y[i] =    11.5500;  t[i++] =   4.0000;
363:   y[i] =    12.0700;  t[i++] =   5.0000;
364:   y[i] =     8.7400;  t[i++] =   6.0000;
365:   y[i] =    80.7000;  t[i++] =    .5000;
366:   y[i] =    61.3000;  t[i++] =    .7500;
367:   y[i] =    47.5000;  t[i++] =   1.0000;
368:    y[i] =   29.0000;  t[i++] =   1.5000;
369:    y[i] =   24.0000;  t[i++] =   2.0000;
370:   y[i] =    17.7000;  t[i++] =   2.5000;
371:   y[i] =    24.5600;  t[i++] =   2.0000;
372:   y[i] =    18.6700;  t[i++] =   2.5000;
373:    y[i] =   16.2400;  t[i++] =   3.0000;
374:   y[i] =     8.7400;  t[i++] =   4.0000;
375:   y[i] =     7.8700;  t[i++] =   5.0000;
376:   y[i] =     8.5100;  t[i++] =   6.0000;
377:   y[i] =    66.7000;  t[i++] =    .5000;
378:   y[i] =    59.2000;  t[i++] =    .7500;
379:   y[i] =    40.8000;  t[i++] =   1.0000;
380:   y[i] =    30.7000;  t[i++] =   1.5000;
381:   y[i] =    25.7000;  t[i++] =   2.0000;
382:   y[i] =    16.3000;  t[i++] =   2.5000;
383:   y[i] =    25.9900;  t[i++] =   2.0000;
384:   y[i] =    16.9500;  t[i++] =   2.5000;
385:   y[i] =    13.3500;  t[i++] =   3.0000;
386:   y[i] =     8.6200;  t[i++] =   4.0000;
387:   y[i] =     7.2000;  t[i++] =   5.0000;
388:   y[i] =     6.6400;  t[i++] =   6.0000;
389:   y[i] =    13.6900;  t[i++] =   3.0000;
390:   y[i] =    81.0000;  t[i++] =    .5000;
391:   y[i] =    64.5000;  t[i++] =    .7500;
392:   y[i] =    35.5000;  t[i++] =   1.5000;
393:    y[i] =   13.3100;  t[i++] =   3.0000;
394:   y[i] =     4.8700;  t[i++] =   6.0000;
395:   y[i] =    12.9400;  t[i++] =   3.0000;
396:   y[i] =     5.0600;  t[i++] =   6.0000;
397:   y[i] =    15.1900;  t[i++] =   3.0000;
398:   y[i] =    14.6200;  t[i++] =   3.0000;
399:   y[i] =    15.6400;  t[i++] =   3.0000;
400:   y[i] =    25.5000;  t[i++] =   1.7500;
401:   y[i] =    25.9500;  t[i++] =   1.7500;
402:   y[i] =    81.7000;  t[i++] =    .5000;
403:   y[i] =    61.6000;  t[i++] =    .7500;
404:   y[i] =    29.8000;  t[i++] =   1.7500;
405:   y[i] =    29.8100;  t[i++] =   1.7500;
406:   y[i] =    17.1700;  t[i++] =   2.7500;
407:   y[i] =    10.3900;  t[i++] =   3.7500;
408:   y[i] =    28.4000;  t[i++] =   1.7500;
409:   y[i] =    28.6900;  t[i++] =   1.7500;
410:   y[i] =    81.3000;  t[i++] =    .5000;
411:   y[i] =    60.9000;  t[i++] =    .7500;
412:   y[i] =    16.6500;  t[i++] =   2.7500;
413:   y[i] =    10.0500;  t[i++] =   3.7500;
414:   y[i] =    28.9000;  t[i++] =   1.7500;
415:   y[i] =    28.9500;  t[i++] =   1.7500;
416:   return(0);
417: }

421: PetscErrorCode MyMonitor(Tao tao, void *ptr)
422: {
423:   PetscReal      fc,gnorm;
424:   PetscInt       its;
425:   PetscViewer    viewer = PETSC_VIEWER_STDOUT_SELF;

429:   TaoGetSolutionStatus(tao,&its,&fc,&gnorm,0,0,0);
430:   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);
431:   ierr=PetscViewerASCIIPrintf(viewer," Function value %g,",(double)fc);
432:   if (gnorm > 1.e-6) {
433:     ierr=PetscViewerASCIIPrintf(viewer," Residual: %g \n",(double)gnorm);
434:   } else if (gnorm > 1.e-11) {
435:     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-6 \n");
436:   } else {
437:     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-11 \n");
438:   }
439:   return(0);
440: }

444: PetscErrorCode TaskWorker(AppCtx *user)
445: {
446:   PetscReal      x[NPARAMETERS],f = 0.0;
447:   PetscMPIInt    tag=IDLE_TAG;
448:   PetscInt       index;
449:   MPI_Status     status;

453:   /* Send check-in message to master */

455:   MPI_Send(&f,1,MPIU_REAL,0,IDLE_TAG,PETSC_COMM_WORLD);
456:   while (tag != DIE_TAG) {
457:     MPI_Recv(x,NPARAMETERS,MPIU_REAL,0,MPI_ANY_TAG,PETSC_COMM_WORLD,&status);
458:     tag = status.MPI_TAG;
459:     if (tag == IDLE_TAG) {
460:       MPI_Send(&f,1,MPIU_REAL,0,IDLE_TAG,PETSC_COMM_WORLD);
461:     } else if (tag != DIE_TAG) {
462:       index = (PetscInt)tag;
463:       ierr=RunSimulation(x,index,&f,user);
464:       ierr=MPI_Send(&f,1,MPIU_REAL,0,tag,PETSC_COMM_WORLD);
465:     }
466:   }
467:   return(0);
468: }

472: PetscErrorCode RunSimulation(PetscReal *x, PetscInt i, PetscReal*f, AppCtx *user)
473: {
474:   PetscReal *t = user->t;
475:   PetscReal *y = user->y;
476:   *f = y[i] - PetscExpScalar(-x[0]*t[i])/(x[1] + x[2]*t[i]);
477:   return(0);
478: }

482: PetscErrorCode StopWorkers(AppCtx *user)
483: {
484:   PetscInt       checkedin;
485:   MPI_Status     status;
486:   PetscReal      f,x[NPARAMETERS];

490:   checkedin=0;
491:   while(checkedin < user->size-1) {
492:     MPI_Recv(&f,1,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PETSC_COMM_WORLD,&status);
493:     checkedin++;
494:     PetscMemzero(x,NPARAMETERS*sizeof(PetscReal));
495:     MPI_Send(x,NPARAMETERS,MPIU_REAL,status.MPI_SOURCE,DIE_TAG,PETSC_COMM_WORLD);
496:   }
497:   return(0);
498: }