1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: posindep.c,v 1.22 1997/12/01 01:56:10 bsmith Exp bsmith $"; 3 #endif 4 /* 5 Code for Timestepping with implicit backwards Euler. 6 */ 7 #include <math.h> 8 #include "src/ts/tsimpl.h" /*I "ts.h" I*/ 9 #include "pinclude/pviewer.h" 10 11 12 typedef struct { 13 Vec update; /* work vector where new solution is formed */ 14 Vec func; /* work vector where F(t[i],u[i]) is stored */ 15 Vec rhs; /* work vector for RHS; vec_sol/dt */ 16 17 /* information used for Pseudo-timestepping */ 18 19 int (*dt)(TS,double*,void*); /* compute next timestep, and related context */ 20 void *dtctx; 21 int (*verify)(TS,Vec,void*,double*,int*); /* verify previous timestep and related context */ 22 void *verifyctx; 23 24 double initial_fnorm,fnorm; /* original and current norm of F(u) */ 25 double fnorm_previous; 26 27 double dt_increment; /* scaling that dt is incremented each time-step */ 28 int increment_dt_from_initial_dt; 29 } TS_Pseudo; 30 31 /* ------------------------------------------------------------------------------*/ 32 #undef __FUNC__ 33 #define __FUNC__ "TSPseudoDefaultTimeStep" 34 /*@C 35 TSPseudoDefaultTimeStep - Default code to compute pseudo-timestepping. 36 Use with TSPseudoSetTimeStep(). 37 38 Input Parameters: 39 . ts - the timestep context 40 . dtctx - unused timestep context 41 42 Output Parameter: 43 . newdt - the timestep to use for the next step 44 45 .keywords: timestep, pseudo, default 46 47 .seealso: TSPseudoSetTimeStep(), TSPseudoComputeTimeStep() 48 @*/ 49 int TSPseudoDefaultTimeStep(TS ts,double* newdt,void* dtctx) 50 { 51 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 52 double inc = pseudo->dt_increment,fnorm_previous = pseudo->fnorm_previous; 53 int ierr; 54 55 PetscFunctionBegin; 56 ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,pseudo->func);CHKERRQ(ierr); 57 ierr = VecNorm(pseudo->func,NORM_2,&pseudo->fnorm); CHKERRQ(ierr); 58 if (pseudo->initial_fnorm == 0.0) { 59 /* first time through so compute initial function norm */ 60 pseudo->initial_fnorm = pseudo->fnorm; 61 fnorm_previous = pseudo->fnorm; 62 } 63 if (pseudo->fnorm == 0.0) { 64 *newdt = 1.e12*inc*ts->time_step; 65 } 66 else if (pseudo->increment_dt_from_initial_dt) { 67 *newdt = inc*ts->initial_time_step*pseudo->initial_fnorm/pseudo->fnorm; 68 } else { 69 *newdt = inc*ts->time_step*fnorm_previous/pseudo->fnorm; 70 } 71 pseudo->fnorm_previous = pseudo->fnorm; 72 PetscFunctionReturn(0); 73 } 74 75 #undef __FUNC__ 76 #define __FUNC__ "TSPseudoSetTimeStep" 77 /*@ 78 TSPseudoSetTimeStep - Sets the user-defined routine to be 79 called at each pseudo-timestep to update the timestep. 80 81 Input Parameters: 82 . ts - timestep context 83 . dt - function to compute timestep 84 . ctx - [optional] user-defined context for private data 85 required by the function (may be PETSC_NULL) 86 87 Calling sequence of func: 88 . func (TS ts,double *newdt,void *ctx); 89 90 . newdt - the newly computed timestep 91 . ctx - [optional] timestep context 92 93 Notes: 94 The routine set here will be called by TSPseudoComputeTimeStep() 95 during the timestepping process. 96 97 .keywords: timestep, pseudo, set 98 99 .seealso: TSPseudoDefaultTimeStep(), TSPseudoComputeTimeStep() 100 @*/ 101 int TSPseudoSetTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 102 { 103 TS_Pseudo *pseudo; 104 105 PetscFunctionBegin; 106 PetscValidHeaderSpecific(ts,TS_COOKIE); 107 if (ts->type != TS_PSEUDO) PetscFunctionReturn(0); 108 109 pseudo = (TS_Pseudo*) ts->data; 110 pseudo->dt = dt; 111 pseudo->dtctx = ctx; 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNC__ 116 #define __FUNC__ "TSPseudoComputeTimeStep" 117 /*@ 118 TSPseudoComputeTimeStep - Computes the next timestep for a currently running 119 pseudo-timestepping process. 120 121 Input Parameter: 122 . ts - timestep context 123 124 Output Parameter: 125 . dt - newly computed timestep 126 127 128 Notes: 129 The routine to be called here to compute the timestep should be 130 set by calling TSPseudoSetTimeStep(). 131 132 .keywords: timestep, pseudo, compute 133 134 .seealso: TSPseudoDefaultTimeStep(), TSPseudoSetTimeStep() 135 @*/ 136 int TSPseudoComputeTimeStep(TS ts,double *dt) 137 { 138 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 139 int ierr; 140 141 PetscFunctionBegin; 142 PLogEventBegin(TS_PseudoComputeTimeStep,ts,0,0,0); 143 ierr = (*pseudo->dt)(ts,dt,pseudo->dtctx); CHKERRQ(ierr); 144 PLogEventEnd(TS_PseudoComputeTimeStep,ts,0,0,0); 145 PetscFunctionReturn(0); 146 } 147 148 149 /* ------------------------------------------------------------------------------*/ 150 #undef __FUNC__ 151 #define __FUNC__ "TSPseudoDefaultVerifyTimeStep" 152 /*@C 153 TSPseudoDefaultVerifyTimeStep - Default code to verify the quality of the last timestep. 154 155 Input Parameters: 156 . ts - the timestep context 157 . dtctx - unused timestep context 158 . update - latest solution vector 159 160 Output Parameters: 161 . newdt - the timestep to use for the next step 162 . flag - flag indicating whether the last time step was acceptable 163 164 Note: 165 This routine always returns a flag of 1, indicating an acceptable 166 timestep. 167 168 .keywords: timestep, pseudo, default, verify 169 170 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoVerifyTimeStep() 171 @*/ 172 int TSPseudoDefaultVerifyTimeStep(TS ts,Vec update,void *dtctx,double *newdt,int *flag) 173 { 174 PetscFunctionBegin; 175 *flag = 1; 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNC__ 180 #define __FUNC__ "TSPseudoSetVerifyTimeStep" 181 /*@ 182 TSPseudoSetVerifyTimeStep - Sets a user-defined routine to verify the quality of the 183 last timestep. 184 185 Input Parameters: 186 . ts - timestep context 187 . dt - user-defined function to verify timestep 188 . ctx - [optional] user-defined context for private data 189 for the timestep verification routine (may be PETSC_NULL) 190 191 Calling sequence of func: 192 . func (TS ts,Vec update,void *ctx,double *newdt,int *flag); 193 194 . update - latest solution vector 195 . ctx - [optional] timestep context 196 . newdt - the timestep to use for the next step 197 . flag - flag indicating whether the last time step was acceptable 198 199 Notes: 200 The routine set here will be called by TSPseudoVerifyTimeStep() 201 during the timestepping process. 202 203 .keywords: timestep, pseudo, set, verify 204 205 .seealso: TSPseudoDefaultVerifyTimeStep(), TSPseudoVerifyTimeStep() 206 @*/ 207 int TSPseudoSetVerifyTimeStep(TS ts,int (*dt)(TS,Vec,void*,double*,int*),void* ctx) 208 { 209 TS_Pseudo *pseudo; 210 211 PetscFunctionBegin; 212 PetscValidHeaderSpecific(ts,TS_COOKIE); 213 if (ts->type != TS_PSEUDO) PetscFunctionReturn(0); 214 215 pseudo = (TS_Pseudo*) ts->data; 216 pseudo->verify = dt; 217 pseudo->verifyctx = ctx; 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNC__ 222 #define __FUNC__ "TSPseudoVerifyTimeStep" 223 /*@ 224 TSPseudoVerifyTimeStep - Verifies whether the last timestep was acceptable. 225 226 Input Parameters: 227 . ts - timestep context 228 . update - latest solution vector 229 230 Output Parameters: 231 . dt - newly computed timestep (if it had to shrink) 232 . flag - indicates if current timestep was ok 233 234 Notes: 235 The routine to be called here to compute the timestep should be 236 set by calling TSPseudoSetVerifyTimeStep(). 237 238 .keywords: timestep, pseudo, verify 239 240 .seealso: TSPseudoSetVerifyTimeStep(), TSPseudoDefaultVerifyTimeStep() 241 @*/ 242 int TSPseudoVerifyTimeStep(TS ts,Vec update,double *dt,int *flag) 243 { 244 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 245 int ierr; 246 247 PetscFunctionBegin; 248 if (!pseudo->verify) {*flag = 1; PetscFunctionReturn(0);} 249 250 ierr = (*pseudo->verify)(ts,update,pseudo->verifyctx,dt,flag ); CHKERRQ(ierr); 251 252 PetscFunctionReturn(0); 253 } 254 255 /* --------------------------------------------------------------------------------*/ 256 257 #undef __FUNC__ 258 #define __FUNC__ "TSStep_Pseudo" 259 static int TSStep_Pseudo(TS ts,int *steps,double *time) 260 { 261 Vec sol = ts->vec_sol; 262 int ierr,i,max_steps = ts->max_steps,its,ok,lits; 263 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 264 double current_time_step; 265 266 PetscFunctionBegin; 267 *steps = -ts->steps; 268 269 ierr = VecCopy(sol,pseudo->update); CHKERRQ(ierr); 270 for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 271 ierr = TSPseudoComputeTimeStep(ts,&ts->time_step); CHKERRQ(ierr); 272 current_time_step = ts->time_step; 273 while (1) { 274 ts->ptime += current_time_step; 275 ierr = SNESSolve(ts->snes,pseudo->update,&its); CHKERRQ(ierr); 276 ierr = SNESGetNumberLinearIterations(ts->snes,&lits); CHKERRQ(ierr); 277 ts->nonlinear_its += PetscAbsInt(its); ts->linear_its += lits; 278 ierr = TSPseudoVerifyTimeStep(ts,pseudo->update,&ts->time_step,&ok); CHKERRQ(ierr); 279 if (ok) break; 280 ts->ptime -= current_time_step; 281 current_time_step = ts->time_step; 282 } 283 ierr = VecCopy(pseudo->update,sol); CHKERRQ(ierr); 284 ts->steps++; 285 ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 286 } 287 288 *steps += ts->steps; 289 *time = ts->ptime; 290 PetscFunctionReturn(0); 291 } 292 293 /*------------------------------------------------------------*/ 294 #undef __FUNC__ 295 #define __FUNC__ "TSDestroy_Pseudo" 296 static int TSDestroy_Pseudo(PetscObject obj ) 297 { 298 TS ts = (TS) obj; 299 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 300 int ierr; 301 302 PetscFunctionBegin; 303 ierr = VecDestroy(pseudo->update); CHKERRQ(ierr); 304 if (pseudo->func) {ierr = VecDestroy(pseudo->func);CHKERRQ(ierr);} 305 if (pseudo->rhs) {ierr = VecDestroy(pseudo->rhs);CHKERRQ(ierr);} 306 if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 307 PetscFree(pseudo); 308 PetscFunctionReturn(0); 309 } 310 311 312 /*------------------------------------------------------------*/ 313 /* 314 This matrix shell multiply where user provided Shell matrix 315 */ 316 317 #undef __FUNC__ 318 #define __FUNC__ "TSPseudoMatMult" 319 int TSPseudoMatMult(Mat mat,Vec x,Vec y) 320 { 321 TS ts; 322 Scalar mdt,mone = -1.0; 323 int ierr; 324 325 PetscFunctionBegin; 326 ierr = MatShellGetContext(mat,(void **)&ts);CHKERRQ(ierr); 327 mdt = 1.0/ts->time_step; 328 329 /* apply user provided function */ 330 ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 331 /* shift and scale by 1/dt - F */ 332 ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 333 PetscFunctionReturn(0); 334 } 335 336 /* 337 This defines the nonlinear equation that is to be solved with SNES 338 339 (U^{n+1} - U^{n})/dt - F(U^{n+1}) 340 */ 341 #undef __FUNC__ 342 #define __FUNC__ "TSPseudoFunction" 343 int TSPseudoFunction(SNES snes,Vec x,Vec y,void *ctx) 344 { 345 TS ts = (TS) ctx; 346 Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 347 int ierr,i,n; 348 349 PetscFunctionBegin; 350 /* apply user provided function */ 351 ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 352 /* compute (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 353 ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 354 ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 355 ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 356 ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 357 for ( i=0; i<n; i++ ) { 358 Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 359 } 360 ierr = VecRestoreArray(ts->vec_sol,&un); 361 ierr = VecRestoreArray(x,&unp1); 362 ierr = VecRestoreArray(y,&Funp1); 363 364 PetscFunctionReturn(0); 365 } 366 367 /* 368 This constructs the Jacobian needed for SNES 369 370 J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 371 */ 372 #undef __FUNC__ 373 #define __FUNC__ "TSPseudoJacobian" 374 int TSPseudoJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 375 { 376 TS ts = (TS) ctx; 377 int ierr; 378 Scalar mone = -1.0, mdt = 1.0/ts->time_step; 379 MatType mtype; 380 381 PetscFunctionBegin; 382 /* construct users Jacobian */ 383 if (ts->rhsjacobian) { 384 ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 385 } 386 387 /* shift and scale Jacobian, if not a shell matrix */ 388 ierr = MatGetType(*AA,&mtype,PETSC_NULL); 389 if (mtype != MATSHELL) { 390 ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 391 ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 392 } 393 ierr = MatGetType(*BB,&mtype,PETSC_NULL); 394 if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 395 ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 396 ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 397 } 398 399 PetscFunctionReturn(0); 400 } 401 402 403 #undef __FUNC__ 404 #define __FUNC__ "TSSetUp_Pseudo" 405 static int TSSetUp_Pseudo(TS ts) 406 { 407 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 408 int ierr, M, m; 409 410 PetscFunctionBegin; 411 ierr = VecDuplicate(ts->vec_sol,&pseudo->update); CHKERRQ(ierr); 412 ierr = VecDuplicate(ts->vec_sol,&pseudo->func); CHKERRQ(ierr); 413 ierr = SNESSetFunction(ts->snes,pseudo->func,TSPseudoFunction,ts);CHKERRQ(ierr); 414 if (ts->Ashell) { /* construct new shell matrix */ 415 ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 416 ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 417 ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 418 ierr = MatShellSetOperation(ts->A,MATOP_MULT,(void*)TSPseudoMatMult);CHKERRQ(ierr); 419 } 420 ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPseudoJacobian,ts);CHKERRQ(ierr); 421 PetscFunctionReturn(0); 422 } 423 /*------------------------------------------------------------*/ 424 425 #undef __FUNC__ 426 #define __FUNC__ "TSPseudoDefaultMonitor" 427 int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 428 { 429 TS_Pseudo *pseudo = (TS_Pseudo*) ts->data; 430 431 PetscFunctionBegin; 432 (*PetscHelpPrintf)(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time,pseudo->fnorm); 433 PetscFunctionReturn(0); 434 } 435 436 #undef __FUNC__ 437 #define __FUNC__ "TSSetFromOptions_Pseudo" 438 static int TSSetFromOptions_Pseudo(TS ts) 439 { 440 int ierr,flg; 441 double inc; 442 443 PetscFunctionBegin; 444 ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 445 446 ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 447 if (flg) { 448 ierr = TSSetMonitor(ts,TSPseudoDefaultMonitor,0); CHKERRQ(ierr); 449 } 450 ierr = OptionsGetDouble(ts->prefix,"-ts_pseudo_increment",&inc,&flg); CHKERRQ(ierr); 451 if (flg) { 452 ierr = TSPseudoSetTimeStepIncrement(ts,inc); CHKERRQ(ierr); 453 } 454 ierr = OptionsHasName(ts->prefix,"-ts_pseudo_increment_dt_from_initial_dt",&flg);CHKERRQ(ierr); 455 if (flg) { 456 ierr = TSPseudoIncrementDtFromInitialDt(ts); CHKERRQ(ierr); 457 } 458 PetscFunctionReturn(0); 459 } 460 461 #undef __FUNC__ 462 #define __FUNC__ "TSPrintHelp_Pseudo" 463 static int TSPrintHelp_Pseudo(TS ts,char *p) 464 { 465 PetscFunctionBegin; 466 (*PetscHelpPrintf)(ts->comm," Options for TS Pseudo timestepper:\n"); 467 (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment <value> : default 1.1\n",p); 468 (*PetscHelpPrintf)(ts->comm," %sts_pseudo_increment_dt_from_initial_dt : use initial_dt *\n",p); 469 (*PetscHelpPrintf)(ts->comm," initial fnorm/current fnorm to determine new timestep\n"); 470 (*PetscHelpPrintf)(ts->comm," default is current_dt * previous fnorm/current fnorm\n"); 471 PetscFunctionReturn(0); 472 } 473 474 #undef __FUNC__ 475 #define __FUNC__ "TSView_Pseudo" 476 static int TSView_Pseudo(PetscObject obj,Viewer viewer) 477 { 478 PetscFunctionBegin; 479 PetscFunctionReturn(0); 480 } 481 482 /* ------------------------------------------------------------ */ 483 #undef __FUNC__ 484 #define __FUNC__ "TSCreate_Pseudo" 485 int TSCreate_Pseudo(TS ts ) 486 { 487 TS_Pseudo *pseudo; 488 int ierr; 489 MatType mtype; 490 491 PetscFunctionBegin; 492 ts->type = TS_PSEUDO; 493 ts->destroy = TSDestroy_Pseudo; 494 ts->printhelp = TSPrintHelp_Pseudo; 495 ts->view = TSView_Pseudo; 496 497 if (ts->problem_type == TS_LINEAR) { 498 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for nonlinear problems"); 499 } 500 if (!ts->A) { 501 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must set Jacobian"); 502 } 503 ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 504 if (mtype == MATSHELL) { 505 ts->Ashell = ts->A; 506 } 507 ts->setup = TSSetUp_Pseudo; 508 ts->step = TSStep_Pseudo; 509 ts->setfromoptions = TSSetFromOptions_Pseudo; 510 511 /* create the required nonlinear solver context */ 512 ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 513 514 pseudo = PetscNew(TS_Pseudo); CHKPTRQ(pseudo); 515 PLogObjectMemory(ts,sizeof(TS_Pseudo)); 516 517 PetscMemzero(pseudo,sizeof(TS_Pseudo)); 518 ts->data = (void *) pseudo; 519 520 pseudo->dt_increment = 1.1; 521 pseudo->increment_dt_from_initial_dt = 0; 522 pseudo->dt = TSPseudoDefaultTimeStep; 523 PetscFunctionReturn(0); 524 } 525 526 527 #undef __FUNC__ 528 #define __FUNC__ "TSPseudoSetTimeStepIncrement" 529 /*@ 530 TSPseudoSetTimeStepIncrement - Sets the scaling increment applied to 531 dt when using the TSPseudoDefaultTimeStep() routine. 532 533 Input Parameters: 534 . ts - the timestep context 535 . inc - the scaling factor >= 1.0 536 537 Options Database Key: 538 $ -ts_pseudo_increment <increment> 539 540 .keywords: timestep, pseudo, set, increment 541 542 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 543 @*/ 544 int TSPseudoSetTimeStepIncrement(TS ts,double inc) 545 { 546 TS_Pseudo *pseudo; 547 548 PetscFunctionBegin; 549 PetscValidHeaderSpecific(ts,TS_COOKIE); 550 if (ts->type != TS_PSEUDO) PetscFunctionReturn(0); 551 552 pseudo = (TS_Pseudo*) ts->data; 553 pseudo->dt_increment = inc; 554 PetscFunctionReturn(0); 555 } 556 557 #undef __FUNC__ 558 #define __FUNC__ "TSPseudoIncrementDtFromInitialDt" 559 /*@ 560 TSPseudoIncrementDtFromInitialDt - Indicates that a new timestep 561 is computed via the formula 562 $ dt = initial_dt*initial_fnorm/current_fnorm 563 rather than the default update, 564 $ dt = current_dt*previous_fnorm/current_fnorm. 565 566 Input Parameter: 567 . ts - the timestep context 568 569 Options Database Key: 570 $ -ts_pseudo_increment_dt_from_initial_dt 571 572 .keywords: timestep, pseudo, set, increment 573 574 .seealso: TSPseudoSetTimeStep(), TSPseudoDefaultTimeStep() 575 @*/ 576 int TSPseudoIncrementDtFromInitialDt(TS ts) 577 { 578 TS_Pseudo *pseudo; 579 580 PetscFunctionBegin; 581 PetscValidHeaderSpecific(ts,TS_COOKIE); 582 if (ts->type != TS_PSEUDO) PetscFunctionReturn(0); 583 584 pseudo = (TS_Pseudo*) ts->data; 585 pseudo->increment_dt_from_initial_dt = 1; 586 PetscFunctionReturn(0); 587 } 588 589 590