1 /* $Id: ts.c,v 1.16 1999/11/10 03:21:40 bsmith Exp bsmith $ */ 2 #include "src/ts/tsimpl.h" /*I "ts.h" I*/ 3 4 #undef __FUNC__ 5 #define __FUNC__ "TSComputeRHSFunction" 6 /* 7 TSComputeRHSFunction - Evaluates the right-hand-side function. 8 9 Note: If the user did not provide a function but merely a matrix, 10 this routine applies the matrix. 11 */ 12 int TSComputeRHSFunction(TS ts,double t,Vec x, Vec y) 13 { 14 int ierr; 15 16 PetscFunctionBegin; 17 PetscValidHeaderSpecific(ts,TS_COOKIE); 18 PetscValidHeader(x); 19 PetscValidHeader(y); 20 21 PLogEventBegin(TS_FunctionEval,ts,x,y,0); 22 if (ts->rhsfunction) { 23 PetscStackPush("TS user right-hand-side function"); 24 ierr = (*ts->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 25 PetscStackPop; 26 } else { 27 if (ts->rhsmatrix) { /* assemble matrix for this timestep */ 28 MatStructure flg; 29 PetscStackPush("TS user right-hand-side matrix function"); 30 ierr = (*ts->rhsmatrix)(ts,t,&ts->A,&ts->B,&flg,ts->jacP);CHKERRQ(ierr); 31 PetscStackPop; 32 } 33 ierr = MatMult(ts->A,x,y);CHKERRQ(ierr); 34 } 35 36 /* apply user-provided boundary conditions (only needed if these are time dependent) */ 37 ierr = TSComputeRHSBoundaryConditions(ts,t,y);CHKERRQ(ierr); 38 PLogEventEnd(TS_FunctionEval,ts,x,y,0); 39 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNC__ 44 #define __FUNC__ "TSSetRHSFunction" 45 /*@C 46 TSSetRHSFunction - Sets the routine for evaluating the function, 47 F(t,u), where U_t = F(t,u). 48 49 Collective on TS 50 51 Input Parameters: 52 + ts - the TS context obtained from TSCreate() 53 . f - routine for evaluating the right-hand-side function 54 - ctx - [optional] user-defined context for private data for the 55 function evaluation routine (may be PETSC_NULL) 56 57 Calling sequence of func: 58 $ func (TS ts,double t,Vec u,Vec F,void *ctx); 59 60 + t - current timestep 61 . u - input vector 62 . F - function vector 63 - ctx - [optional] user-defined function context 64 65 Important: 66 The user MUST call either this routine or TSSetRHSMatrix(). 67 68 Level: beginner 69 70 .keywords: TS, timestep, set, right-hand-side, function 71 72 .seealso: TSSetRHSMatrix() 73 @*/ 74 int TSSetRHSFunction(TS ts,int (*f)(TS,double,Vec,Vec,void*),void *ctx) 75 { 76 PetscFunctionBegin; 77 78 PetscValidHeaderSpecific(ts,TS_COOKIE); 79 if (ts->problem_type == TS_LINEAR) { 80 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Cannot set function for linear problem"); 81 } 82 ts->rhsfunction = f; 83 ts->funP = ctx; 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNC__ 88 #define __FUNC__ "TSSetRHSMatrix" 89 /*@C 90 TSSetRHSMatrix - Sets the function to compute the matrix A, where U_t = A(t) U. 91 Also sets the location to store A. 92 93 Collective on TS 94 95 Input Parameters: 96 + ts - the TS context obtained from TSCreate() 97 . A - matrix 98 . B - preconditioner matrix (usually same as A) 99 . f - the matrix evaluation routine; use PETSC_NULL (PETSC_NULL_FUNCTION in fortran) 100 if A is not a function of t. 101 - ctx - [optional] user-defined context for private data for the 102 matrix evaluation routine (may be PETSC_NULL) 103 104 Calling sequence of func: 105 $ func (TS ts,double t,Mat *A,Mat *B,int *flag,void *ctx); 106 107 + t - current timestep 108 . A - matrix A, where U_t = A(t) U 109 . B - preconditioner matrix, usually the same as A 110 . flag - flag indicating information about the preconditioner matrix 111 structure (same as flag in SLESSetOperators()) 112 - ctx - [optional] user-defined context for matrix evaluation routine 113 114 Notes: 115 See SLESSetOperators() for important information about setting the flag 116 output parameter in the routine func(). Be sure to read this information! 117 118 The routine func() takes Mat * as the matrix arguments rather than Mat. 119 This allows the matrix evaluation routine to replace A and/or B with a 120 completely new new matrix structure (not just different matrix elements) 121 when appropriate, for instance, if the nonzero structure is changing 122 throughout the global iterations. 123 124 Important: 125 The user MUST call either this routine or TSSetRHSFunction(). 126 127 Level: beginner 128 129 .keywords: TS, timestep, set, right-hand-side, matrix 130 131 .seealso: TSSetRHSFunction() 132 @*/ 133 int TSSetRHSMatrix(TS ts,Mat A, Mat B,int (*f)(TS,double,Mat*,Mat*,MatStructure*,void*),void *ctx) 134 { 135 PetscFunctionBegin; 136 PetscValidHeaderSpecific(ts,TS_COOKIE); 137 PetscValidHeaderSpecific(A,MAT_COOKIE); 138 PetscValidHeaderSpecific(B,MAT_COOKIE); 139 PetscCheckSameComm(ts,A); 140 PetscCheckSameComm(ts,B); 141 if (ts->problem_type == TS_NONLINEAR) { 142 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for nonlinear problems; use TSSetRHSJacobian()"); 143 } 144 145 ts->rhsmatrix = f; 146 ts->jacP = ctx; 147 ts->A = A; 148 ts->B = B; 149 150 PetscFunctionReturn(0); 151 } 152 153 #undef __FUNC__ 154 #define __FUNC__ "TSSetRHSJacobian" 155 /*@C 156 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 157 where U_t = F(U,t), as well as the location to store the matrix. 158 159 Collective on TS 160 161 Input Parameters: 162 + ts - the TS context obtained from TSCreate() 163 . A - Jacobian matrix 164 . B - preconditioner matrix (usually same as A) 165 . f - the Jacobian evaluation routine 166 - ctx - [optional] user-defined context for private data for the 167 Jacobian evaluation routine (may be PETSC_NULL) 168 169 Calling sequence of func: 170 $ func (TS ts,double t,Vec u,Mat *A,Mat *B,int *flag,void *ctx); 171 172 + t - current timestep 173 . u - input vector 174 . A - matrix A, where U_t = A(t)u 175 . B - preconditioner matrix, usually the same as A 176 . flag - flag indicating information about the preconditioner matrix 177 structure (same as flag in SLESSetOperators()) 178 - ctx - [optional] user-defined context for matrix evaluation routine 179 180 Notes: 181 See SLESSetOperators() for important information about setting the flag 182 output parameter in the routine func(). Be sure to read this information! 183 184 The routine func() takes Mat * as the matrix arguments rather than Mat. 185 This allows the matrix evaluation routine to replace A and/or B with a 186 completely new new matrix structure (not just different matrix elements) 187 when appropriate, for instance, if the nonzero structure is changing 188 throughout the global iterations. 189 190 Level: beginner 191 192 .keywords: TS, timestep, set, right-hand-side, Jacobian 193 194 .seealso: TSDefaultComputeJacobianColor(), 195 SNESDefaultComputeJacobianColor() 196 197 @*/ 198 int TSSetRHSJacobian(TS ts,Mat A, Mat B,int (*f)(TS,double,Vec,Mat*,Mat*, 199 MatStructure*,void*),void *ctx) 200 { 201 PetscFunctionBegin; 202 PetscValidHeaderSpecific(ts,TS_COOKIE); 203 PetscValidHeaderSpecific(A,MAT_COOKIE); 204 PetscValidHeaderSpecific(B,MAT_COOKIE); 205 PetscCheckSameComm(ts,A); 206 PetscCheckSameComm(ts,B); 207 if (ts->problem_type != TS_NONLINEAR) { 208 SETERRQ(PETSC_ERR_ARG_WRONG,0,"Not for linear problems; use TSSetRHSMatrix()"); 209 } 210 211 ts->rhsjacobian = f; 212 ts->jacP = ctx; 213 ts->A = A; 214 ts->B = B; 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNC__ 219 #define __FUNC__ "TSComputeRHSBoundaryConditions" 220 /* 221 TSComputeRHSBoundaryConditions - Evaluates the boundary condition function. 222 223 Note: If the user did not provide a function but merely a matrix, 224 this routine applies the matrix. 225 */ 226 int TSComputeRHSBoundaryConditions(TS ts,double t,Vec x) 227 { 228 int ierr; 229 230 PetscFunctionBegin; 231 PetscValidHeaderSpecific(ts,TS_COOKIE); 232 PetscValidHeader(x); 233 PetscCheckSameComm(ts,x); 234 235 if (ts->rhsbc) { 236 PetscStackPush("TS user boundary condition function"); 237 ierr = (*ts->rhsbc)(ts,t,x,ts->bcP);CHKERRQ(ierr); 238 PetscStackPop; 239 PetscFunctionReturn(0); 240 } 241 242 PetscFunctionReturn(0); 243 } 244 245 #undef __FUNC__ 246 #define __FUNC__ "TSSetRHSBoundaryConditions" 247 /*@C 248 TSSetRHSBoundaryConditions - Sets the routine for evaluating the function, 249 boundary conditions for the function F. 250 251 Collective on TS 252 253 Input Parameters: 254 + ts - the TS context obtained from TSCreate() 255 . f - routine for evaluating the boundary condition function 256 - ctx - [optional] user-defined context for private data for the 257 function evaluation routine (may be PETSC_NULL) 258 259 Calling sequence of func: 260 $ func (TS ts,double t,Vec F,void *ctx); 261 262 + t - current timestep 263 . F - function vector 264 - ctx - [optional] user-defined function context 265 266 Level: intermediate 267 268 .keywords: TS, timestep, set, boundary conditions, function 269 @*/ 270 int TSSetRHSBoundaryConditions(TS ts,int (*f)(TS,double,Vec,void*),void *ctx) 271 { 272 PetscFunctionBegin; 273 274 PetscValidHeaderSpecific(ts,TS_COOKIE); 275 if (ts->problem_type != TS_LINEAR) { 276 SETERRQ(PETSC_ERR_ARG_WRONG,0,"For linear problems only"); 277 } 278 ts->rhsbc = f; 279 ts->bcP = ctx; 280 PetscFunctionReturn(0); 281 } 282 283 #undef __FUNC__ 284 #define __FUNC__ "TSView" 285 /*@ 286 TSView - Prints the TS data structure. 287 288 Collective on TS, unless Viewer is VIEWER_STDOUT_SELF 289 290 Input Parameters: 291 + ts - the TS context obtained from TSCreate() 292 - viewer - visualization context 293 294 Options Database Key: 295 . -ts_view - calls TSView() at end of TSStep() 296 297 Notes: 298 The available visualization contexts include 299 + VIEWER_STDOUT_SELF - standard output (default) 300 - VIEWER_STDOUT_WORLD - synchronized standard 301 output where only the first processor opens 302 the file. All other processors send their 303 data to the first processor to print. 304 305 The user can open an alternative visualization context with 306 ViewerASCIIOpen() - output to a specified file. 307 308 Level: beginner 309 310 .keywords: TS, timestep, view 311 312 .seealso: ViewerASCIIOpen() 313 @*/ 314 int TSView(TS ts,Viewer viewer) 315 { 316 int ierr; 317 char *type; 318 PetscTruth isascii,isstring; 319 320 PetscFunctionBegin; 321 PetscValidHeaderSpecific(ts,TS_COOKIE); 322 if (!viewer) viewer = VIEWER_STDOUT_SELF; 323 PetscValidHeaderSpecific(viewer,VIEWER_COOKIE); 324 PetscCheckSameComm(ts,viewer); 325 326 ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 327 ierr = PetscTypeCompare((PetscObject)viewer,STRING_VIEWER,&isstring);CHKERRQ(ierr); 328 if (isascii) { 329 ierr = ViewerASCIIPrintf(viewer,"TS Object:\n");CHKERRQ(ierr); 330 ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr); 331 if (type) { 332 ierr = ViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 333 } else { 334 ierr = ViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 335 } 336 if (ts->view) { 337 ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 338 ierr = (*ts->view)(ts,viewer);CHKERRQ(ierr); 339 ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 340 } 341 ierr = ViewerASCIIPrintf(viewer," maximum steps=%d\n",ts->max_steps);CHKERRQ(ierr); 342 ierr = ViewerASCIIPrintf(viewer," maximum time=%g\n",ts->max_time);CHKERRQ(ierr); 343 if (ts->problem_type == TS_NONLINEAR) { 344 ierr = ViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%d\n",ts->nonlinear_its);CHKERRQ(ierr); 345 } 346 ierr = ViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",ts->linear_its);CHKERRQ(ierr); 347 } else if (isstring) { 348 ierr = TSGetType(ts,(TSType *)&type);CHKERRQ(ierr); 349 ierr = ViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 350 } 351 ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 352 if (ts->sles) {ierr = SLESView(ts->sles,viewer);CHKERRQ(ierr);} 353 if (ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 354 ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 355 PetscFunctionReturn(0); 356 } 357 358 359 #undef __FUNC__ 360 #define __FUNC__ "TSSetApplicationContext" 361 /*@C 362 TSSetApplicationContext - Sets an optional user-defined context for 363 the timesteppers. 364 365 Collective on TS 366 367 Input Parameters: 368 + ts - the TS context obtained from TSCreate() 369 - usrP - optional user context 370 371 Level: intermediate 372 373 .keywords: TS, timestep, set, application, context 374 375 .seealso: TSGetApplicationContext() 376 @*/ 377 int TSSetApplicationContext(TS ts,void *usrP) 378 { 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(ts,TS_COOKIE); 381 ts->user = usrP; 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNC__ 386 #define __FUNC__ "TSGetApplicationContext" 387 /*@C 388 TSGetApplicationContext - Gets the user-defined context for the 389 timestepper. 390 391 Not Collective 392 393 Input Parameter: 394 . ts - the TS context obtained from TSCreate() 395 396 Output Parameter: 397 . usrP - user context 398 399 Level: intermediate 400 401 .keywords: TS, timestep, get, application, context 402 403 .seealso: TSSetApplicationContext() 404 @*/ 405 int TSGetApplicationContext( TS ts, void **usrP ) 406 { 407 PetscFunctionBegin; 408 PetscValidHeaderSpecific(ts,TS_COOKIE); 409 *usrP = ts->user; 410 PetscFunctionReturn(0); 411 } 412 413 #undef __FUNC__ 414 #define __FUNC__ "TSGetTimeStepNumber" 415 /*@ 416 TSGetTimeStepNumber - Gets the current number of timesteps. 417 418 Not Collective 419 420 Input Parameter: 421 . ts - the TS context obtained from TSCreate() 422 423 Output Parameter: 424 . iter - number steps so far 425 426 Level: intermediate 427 428 .keywords: TS, timestep, get, iteration, number 429 @*/ 430 int TSGetTimeStepNumber(TS ts,int* iter) 431 { 432 PetscFunctionBegin; 433 PetscValidHeaderSpecific(ts,TS_COOKIE); 434 *iter = ts->steps; 435 PetscFunctionReturn(0); 436 } 437 438 #undef __FUNC__ 439 #define __FUNC__ "TSSetInitialTimeStep" 440 /*@ 441 TSSetInitialTimeStep - Sets the initial timestep to be used, 442 as well as the initial time. 443 444 Collective on TS 445 446 Input Parameters: 447 + ts - the TS context obtained from TSCreate() 448 . initial_time - the initial time 449 - time_step - the size of the timestep 450 451 Level: intermediate 452 453 .seealso: TSSetTimeStep(), TSGetTimeStep() 454 455 .keywords: TS, set, initial, timestep 456 @*/ 457 int TSSetInitialTimeStep(TS ts,double initial_time,double time_step) 458 { 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(ts,TS_COOKIE); 461 ts->time_step = time_step; 462 ts->initial_time_step = time_step; 463 ts->ptime = initial_time; 464 PetscFunctionReturn(0); 465 } 466 467 #undef __FUNC__ 468 #define __FUNC__ "TSSetTimeStep" 469 /*@ 470 TSSetTimeStep - Allows one to reset the timestep at any time, 471 useful for simple pseudo-timestepping codes. 472 473 Collective on TS 474 475 Input Parameters: 476 + ts - the TS context obtained from TSCreate() 477 - time_step - the size of the timestep 478 479 Level: intermediate 480 481 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 482 483 .keywords: TS, set, timestep 484 @*/ 485 int TSSetTimeStep(TS ts,double time_step) 486 { 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(ts,TS_COOKIE); 489 ts->time_step = time_step; 490 PetscFunctionReturn(0); 491 } 492 493 #undef __FUNC__ 494 #define __FUNC__ "TSGetTimeStep" 495 /*@ 496 TSGetTimeStep - Gets the current timestep size. 497 498 Not Collective 499 500 Input Parameter: 501 . ts - the TS context obtained from TSCreate() 502 503 Output Parameter: 504 . dt - the current timestep size 505 506 Level: intermediate 507 508 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 509 510 .keywords: TS, get, timestep 511 @*/ 512 int TSGetTimeStep(TS ts,double* dt) 513 { 514 PetscFunctionBegin; 515 PetscValidHeaderSpecific(ts,TS_COOKIE); 516 *dt = ts->time_step; 517 PetscFunctionReturn(0); 518 } 519 520 #undef __FUNC__ 521 #define __FUNC__ "TSGetSolution" 522 /*@C 523 TSGetSolution - Returns the solution at the present timestep. It 524 is valid to call this routine inside the function that you are evaluating 525 in order to move to the new timestep. This vector not changed until 526 the solution at the next timestep has been calculated. 527 528 Not Collective, but Vec returned is parallel if TS is parallel 529 530 Input Parameter: 531 . ts - the TS context obtained from TSCreate() 532 533 Output Parameter: 534 . v - the vector containing the solution 535 536 Level: intermediate 537 538 .seealso: TSGetTimeStep() 539 540 .keywords: TS, timestep, get, solution 541 @*/ 542 int TSGetSolution(TS ts,Vec *v) 543 { 544 PetscFunctionBegin; 545 PetscValidHeaderSpecific(ts,TS_COOKIE); 546 *v = ts->vec_sol_always; 547 PetscFunctionReturn(0); 548 } 549 550 #undef __FUNC__ 551 #define __FUNC__ "TSPublish_Petsc" 552 static int TSPublish_Petsc(PetscObject obj) 553 { 554 #if defined(PETSC_HAVE_AMS) 555 TS v = (TS) obj; 556 int ierr; 557 #endif 558 559 PetscFunctionBegin; 560 561 #if defined(PETSC_HAVE_AMS) 562 /* if it is already published then return */ 563 if (v->amem >=0 ) PetscFunctionReturn(0); 564 565 ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 566 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Step",&v->steps,1,AMS_INT,AMS_READ, 567 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 568 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Time",&v->ptime,1,AMS_DOUBLE,AMS_READ, 569 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 570 ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"CurrentTimeStep",&v->time_step,1, 571 AMS_DOUBLE,AMS_READ,AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 572 ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 573 #endif 574 PetscFunctionReturn(0); 575 } 576 577 /* -----------------------------------------------------------*/ 578 579 #undef __FUNC__ 580 #define __FUNC__ "TSCreate" 581 /*@C 582 TSCreate - Creates a timestepper context. 583 584 Collective on MPI_Comm 585 586 Input Parameter: 587 + comm - MPI communicator 588 - type - One of TS_LINEAR,TS_NONLINEAR 589 where these types refer to problems of the forms 590 .vb 591 U_t = A U 592 U_t = A(t) U 593 U_t = F(t,U) 594 .ve 595 596 Output Parameter: 597 . outts - the new TS context 598 599 Level: beginner 600 601 .keywords: TS, timestep, create, context 602 603 .seealso: TSSetUp(), TSStep(), TSDestroy() 604 @*/ 605 int TSCreate(MPI_Comm comm,TSProblemType problemtype,TS *outts) 606 { 607 TS ts; 608 609 PetscFunctionBegin; 610 *outts = 0; 611 PetscHeaderCreate(ts,_p_TS,int,TS_COOKIE,-1,"TS",comm,TSDestroy,TSView); 612 PLogObjectCreate(ts); 613 ts->bops->publish = TSPublish_Petsc; 614 ts->max_steps = 5000; 615 ts->max_time = 5.0; 616 ts->time_step = .1; 617 ts->initial_time_step = ts->time_step; 618 ts->steps = 0; 619 ts->ptime = 0.0; 620 ts->data = 0; 621 ts->view = 0; 622 ts->setupcalled = 0; 623 ts->problem_type = problemtype; 624 ts->numbermonitors = 0; 625 ts->linear_its = 0; 626 ts->nonlinear_its = 0; 627 628 *outts = ts; 629 PetscFunctionReturn(0); 630 } 631 632 /* ----- Routines to initialize and destroy a timestepper ---- */ 633 634 #undef __FUNC__ 635 #define __FUNC__ "TSSetUp" 636 /*@ 637 TSSetUp - Sets up the internal data structures for the later use 638 of a timestepper. 639 640 Collective on TS 641 642 Input Parameter: 643 . ts - the TS context obtained from TSCreate() 644 645 Notes: 646 For basic use of the TS solvers the user need not explicitly call 647 TSSetUp(), since these actions will automatically occur during 648 the call to TSStep(). However, if one wishes to control this 649 phase separately, TSSetUp() should be called after TSCreate() 650 and optional routines of the form TSSetXXX(), but before TSStep(). 651 652 Level: advanced 653 654 .keywords: TS, timestep, setup 655 656 .seealso: TSCreate(), TSStep(), TSDestroy() 657 @*/ 658 int TSSetUp(TS ts) 659 { 660 int ierr; 661 662 PetscFunctionBegin; 663 PetscValidHeaderSpecific(ts,TS_COOKIE); 664 if (!ts->vec_sol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call TSSetSolution() first"); 665 if (!ts->type_name) { 666 ierr = TSSetType(ts,TS_EULER);CHKERRQ(ierr); 667 } 668 ierr = (*ts->setup)(ts);CHKERRQ(ierr); 669 ts->setupcalled = 1; 670 PetscFunctionReturn(0); 671 } 672 673 #undef __FUNC__ 674 #define __FUNC__ "TSDestroy" 675 /*@C 676 TSDestroy - Destroys the timestepper context that was created 677 with TSCreate(). 678 679 Collective on TS 680 681 Input Parameter: 682 . ts - the TS context obtained from TSCreate() 683 684 Level: beginner 685 686 .keywords: TS, timestepper, destroy 687 688 .seealso: TSCreate(), TSSetUp(), TSSolve() 689 @*/ 690 int TSDestroy(TS ts) 691 { 692 int ierr; 693 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(ts,TS_COOKIE); 696 if (--ts->refct > 0) PetscFunctionReturn(0); 697 698 /* if memory was published with AMS then destroy it */ 699 ierr = PetscObjectDepublish(ts);CHKERRQ(ierr); 700 701 if (ts->sles) {ierr = SLESDestroy(ts->sles);CHKERRQ(ierr);} 702 if (ts->snes) {ierr = SNESDestroy(ts->snes);CHKERRQ(ierr);} 703 ierr = (*(ts)->destroy)(ts);CHKERRQ(ierr); 704 PLogObjectDestroy((PetscObject)ts); 705 PetscHeaderDestroy((PetscObject)ts); 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNC__ 710 #define __FUNC__ "TSGetSNES" 711 /*@C 712 TSGetSNES - Returns the SNES (nonlinear solver) associated with 713 a TS (timestepper) context. Valid only for nonlinear problems. 714 715 Not Collective, but SNES is parallel if TS is parallel 716 717 Input Parameter: 718 . ts - the TS context obtained from TSCreate() 719 720 Output Parameter: 721 . snes - the nonlinear solver context 722 723 Notes: 724 The user can then directly manipulate the SNES context to set various 725 options, etc. Likewise, the user can then extract and manipulate the 726 SLES, KSP, and PC contexts as well. 727 728 TSGetSNES() does not work for integrators that do not use SNES; in 729 this case TSGetSNES() returns PETSC_NULL in snes. 730 731 Level: beginner 732 733 .keywords: timestep, get, SNES 734 @*/ 735 int TSGetSNES(TS ts,SNES *snes) 736 { 737 PetscFunctionBegin; 738 PetscValidHeaderSpecific(ts,TS_COOKIE); 739 if (ts->problem_type == TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Nonlinear only; use TSGetSLES()"); 740 *snes = ts->snes; 741 PetscFunctionReturn(0); 742 } 743 744 #undef __FUNC__ 745 #define __FUNC__ "TSGetSLES" 746 /*@C 747 TSGetSLES - Returns the SLES (linear solver) associated with 748 a TS (timestepper) context. 749 750 Not Collective, but SLES is parallel if TS is parallel 751 752 Input Parameter: 753 . ts - the TS context obtained from TSCreate() 754 755 Output Parameter: 756 . sles - the nonlinear solver context 757 758 Notes: 759 The user can then directly manipulate the SLES context to set various 760 options, etc. Likewise, the user can then extract and manipulate the 761 KSP and PC contexts as well. 762 763 TSGetSLES() does not work for integrators that do not use SLES; 764 in this case TSGetSLES() returns PETSC_NULL in sles. 765 766 Level: beginner 767 768 .keywords: timestep, get, SLES 769 @*/ 770 int TSGetSLES(TS ts,SLES *sles) 771 { 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(ts,TS_COOKIE); 774 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Linear only; use TSGetSNES()"); 775 *sles = ts->sles; 776 PetscFunctionReturn(0); 777 } 778 779 /* ----------- Routines to set solver parameters ---------- */ 780 781 #undef __FUNC__ 782 #define __FUNC__ "TSSetDuration" 783 /*@ 784 TSSetDuration - Sets the maximum number of timesteps to use and 785 maximum time for iteration. 786 787 Collective on TS 788 789 Input Parameters: 790 + ts - the TS context obtained from TSCreate() 791 . maxsteps - maximum number of iterations to use 792 - maxtime - final time to iterate to 793 794 Options Database Keys: 795 . -ts_max_steps <maxsteps> - Sets maxsteps 796 . -ts_max_time <maxtime> - Sets maxtime 797 798 Notes: 799 The default maximum number of iterations is 5000. Default time is 5.0 800 801 Level: intermediate 802 803 .keywords: TS, timestep, set, maximum, iterations 804 @*/ 805 int TSSetDuration(TS ts,int maxsteps,double maxtime) 806 { 807 PetscFunctionBegin; 808 PetscValidHeaderSpecific(ts,TS_COOKIE); 809 ts->max_steps = maxsteps; 810 ts->max_time = maxtime; 811 PetscFunctionReturn(0); 812 } 813 814 #undef __FUNC__ 815 #define __FUNC__ "TSSetSolution" 816 /*@ 817 TSSetSolution - Sets the initial solution vector 818 for use by the TS routines. 819 820 Collective on TS and Vec 821 822 Input Parameters: 823 + ts - the TS context obtained from TSCreate() 824 - x - the solution vector 825 826 Level: beginner 827 828 .keywords: TS, timestep, set, solution, initial conditions 829 @*/ 830 int TSSetSolution(TS ts,Vec x) 831 { 832 PetscFunctionBegin; 833 PetscValidHeaderSpecific(ts,TS_COOKIE); 834 ts->vec_sol = ts->vec_sol_always = x; 835 PetscFunctionReturn(0); 836 } 837 838 /* ------------ Routines to set performance monitoring options ----------- */ 839 840 #undef __FUNC__ 841 #define __FUNC__ "TSSetMonitor" 842 /*@C 843 TSSetMonitor - Sets an ADDITIONAL function that is to be used at every 844 timestep to display the iteration's progress. 845 846 Collective on TS 847 848 Input Parameters: 849 + ts - the TS context obtained from TSCreate() 850 . func - monitoring routine 851 - mctx - [optional] user-defined context for private data for the 852 monitor routine (may be PETSC_NULL) 853 854 Calling sequence of func: 855 $ int func(TS ts,int steps,double time,Vec x,void *mctx) 856 857 + ts - the TS context 858 . steps - iteration number 859 . time - current timestep 860 . x - current iterate 861 - mctx - [optional] monitoring context 862 863 Notes: 864 This routine adds an additional monitor to the list of monitors that 865 already has been loaded. 866 867 Level: intermediate 868 869 .keywords: TS, timestep, set, monitor 870 871 .seealso: TSDefaultMonitor(), TSClearMonitor() 872 @*/ 873 int TSSetMonitor(TS ts, int (*monitor)(TS,int,double,Vec,void*), void *mctx ) 874 { 875 PetscFunctionBegin; 876 PetscValidHeaderSpecific(ts,TS_COOKIE); 877 if (ts->numbermonitors >= MAXTSMONITORS) { 878 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 879 } 880 ts->monitor[ts->numbermonitors] = monitor; 881 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 882 PetscFunctionReturn(0); 883 } 884 885 #undef __FUNC__ 886 #define __FUNC__ "TSClearMonitor" 887 /*@C 888 TSClearMonitor - Clears all the monitors that have been set on a time-step object. 889 890 Collective on TS 891 892 Input Parameters: 893 . ts - the TS context obtained from TSCreate() 894 895 Notes: 896 There is no way to remove a single, specific monitor. 897 898 Level: intermediate 899 900 .keywords: TS, timestep, set, monitor 901 902 .seealso: TSDefaultMonitor(), TSSetMonitor() 903 @*/ 904 int TSClearMonitor(TS ts) 905 { 906 PetscFunctionBegin; 907 PetscValidHeaderSpecific(ts,TS_COOKIE); 908 ts->numbermonitors = 0; 909 PetscFunctionReturn(0); 910 } 911 912 #undef __FUNC__ 913 #define __FUNC__ "TSDefaultMonitor" 914 int TSDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 915 { 916 int ierr; 917 918 PetscFunctionBegin; 919 ierr = PetscPrintf(ts->comm,"timestep %d dt %g time %g\n",step,ts->time_step,time);CHKERRQ(ierr); 920 PetscFunctionReturn(0); 921 } 922 923 #undef __FUNC__ 924 #define __FUNC__ "TSStep" 925 /*@ 926 TSStep - Steps the requested number of timesteps. 927 928 Collective on TS 929 930 Input Parameter: 931 . ts - the TS context obtained from TSCreate() 932 933 Output Parameters: 934 + steps - number of iterations until termination 935 - time - time until termination 936 937 Level: beginner 938 939 .keywords: TS, timestep, solve 940 941 .seealso: TSCreate(), TSSetUp(), TSDestroy() 942 @*/ 943 int TSStep(TS ts,int *steps,double *time) 944 { 945 int ierr; 946 PetscTruth flg; 947 948 PetscFunctionBegin; 949 PetscValidHeaderSpecific(ts,TS_COOKIE); 950 if (!ts->setupcalled) {ierr = TSSetUp(ts);CHKERRQ(ierr);} 951 PLogEventBegin(TS_Step,ts,0,0,0); 952 ierr = (*(ts)->step)(ts,steps,time);CHKERRQ(ierr); 953 PLogEventEnd(TS_Step,ts,0,0,0); 954 ierr = OptionsHasName(ts->prefix,"-ts_view",&flg);CHKERRQ(ierr); 955 if (flg) { 956 ierr = TSView(ts,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 957 } 958 PetscFunctionReturn(0); 959 } 960 961 #undef __FUNC__ 962 #define __FUNC__ "TSMonitor" 963 /* 964 Runs the user provided monitor routines, if they exists. 965 */ 966 int TSMonitor(TS ts,int step,double time,Vec x) 967 { 968 int i,ierr,n = ts->numbermonitors; 969 970 PetscFunctionBegin; 971 for ( i=0; i<n; i++ ) { 972 ierr = (*ts->monitor[i])(ts,step,time,x,ts->monitorcontext[i]);CHKERRQ(ierr); 973 } 974 PetscFunctionReturn(0); 975 } 976 977 /* ------------------------------------------------------------------------*/ 978 979 #undef __FUNC__ 980 #define __FUNC__ "TSLGMonitorCreate" 981 /*@C 982 TSLGMonitorCreate - Creates a line graph context for use with 983 TS to monitor convergence of preconditioned residual norms. 984 985 Collective on TS 986 987 Input Parameters: 988 + host - the X display to open, or null for the local machine 989 . label - the title to put in the title bar 990 . x, y - the screen coordinates of the upper left coordinate of the window 991 - m, n - the screen width and height in pixels 992 993 Output Parameter: 994 . draw - the drawing context 995 996 Options Database Key: 997 . -ts_xmonitor - automatically sets line graph monitor 998 999 Notes: 1000 Use TSLGMonitorDestroy() to destroy this line graph, not DrawLGDestroy(). 1001 1002 Level: intermediate 1003 1004 .keywords: TS, monitor, line graph, residual, seealso 1005 1006 .seealso: TSLGMonitorDestroy(), TSSetMonitor() 1007 1008 @*/ 1009 int TSLGMonitorCreate(char *host,char *label,int x,int y,int m,int n, DrawLG *draw) 1010 { 1011 Draw win; 1012 int ierr; 1013 1014 PetscFunctionBegin; 1015 ierr = DrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1016 ierr = DrawSetType(win,DRAW_X);CHKERRQ(ierr); 1017 ierr = DrawLGCreate(win,1,draw);CHKERRQ(ierr); 1018 ierr = DrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1019 1020 PLogObjectParent(*draw,win); 1021 PetscFunctionReturn(0); 1022 } 1023 1024 #undef __FUNC__ 1025 #define __FUNC__ "TSLGMonitor" 1026 int TSLGMonitor(TS ts,int n,double time,Vec v,void *monctx) 1027 { 1028 DrawLG lg = (DrawLG) monctx; 1029 double x,y = time; 1030 int ierr; 1031 1032 PetscFunctionBegin; 1033 if (!monctx) { 1034 MPI_Comm comm; 1035 Viewer viewer; 1036 1037 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1038 viewer = VIEWER_DRAW_(comm); 1039 ierr = ViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1040 } 1041 1042 if (!n) {ierr = DrawLGReset(lg);CHKERRQ(ierr);} 1043 x = (double) n; 1044 ierr = DrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1045 if (n < 20 || (n % 5)) { 1046 ierr = DrawLGDraw(lg);CHKERRQ(ierr); 1047 } 1048 PetscFunctionReturn(0); 1049 } 1050 1051 #undef __FUNC__ 1052 #define __FUNC__ "TSLGMonitorDestroy" 1053 /*@C 1054 TSLGMonitorDestroy - Destroys a line graph context that was created 1055 with TSLGMonitorCreate(). 1056 1057 Collective on DrawLG 1058 1059 Input Parameter: 1060 . draw - the drawing context 1061 1062 Level: intermediate 1063 1064 .keywords: TS, monitor, line graph, destroy 1065 1066 .seealso: TSLGMonitorCreate(), TSSetMonitor(), TSLGMonitor(); 1067 @*/ 1068 int TSLGMonitorDestroy(DrawLG drawlg) 1069 { 1070 Draw draw; 1071 int ierr; 1072 1073 PetscFunctionBegin; 1074 ierr = DrawLGGetDraw(drawlg,&draw);CHKERRQ(ierr); 1075 ierr = DrawDestroy(draw);CHKERRQ(ierr); 1076 ierr = DrawLGDestroy(drawlg);CHKERRQ(ierr); 1077 PetscFunctionReturn(0); 1078 } 1079 1080 #undef __FUNC__ 1081 #define __FUNC__ "TSGetTime" 1082 /*@ 1083 TSGetTime - Gets the current time. 1084 1085 Not Collective 1086 1087 Input Parameter: 1088 . ts - the TS context obtained from TSCreate() 1089 1090 Output Parameter: 1091 . t - the current time 1092 1093 Contributed by: Matthew Knepley 1094 1095 Level: beginner 1096 1097 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1098 1099 .keywords: TS, get, time 1100 @*/ 1101 int TSGetTime(TS ts, double* t) 1102 { 1103 PetscFunctionBegin; 1104 PetscValidHeaderSpecific(ts, TS_COOKIE); 1105 *t = ts->ptime; 1106 PetscFunctionReturn(0); 1107 } 1108 1109 #undef __FUNC__ 1110 #define __FUNC__ "TSGetProblemType" 1111 /*@C 1112 TSGetProblemType - Returns the problem type of a TS (timestepper) context. 1113 1114 Not Collective 1115 1116 Input Parameter: 1117 . ts - The TS context obtained from TSCreate() 1118 1119 Output Parameter: 1120 . type - The problem type, TS_LINEAR or TS_NONLINEAR 1121 1122 Level: intermediate 1123 1124 Contributed by: Matthew Knepley 1125 1126 .keywords: ts, get, type 1127 1128 @*/ 1129 int TSGetProblemType(TS ts, TSProblemType *type) 1130 { 1131 PetscFunctionBegin; 1132 PetscValidHeaderSpecific(ts, TS_COOKIE); 1133 *type = ts->problem_type; 1134 PetscFunctionReturn(0); 1135 } 1136 1137 #undef __FUNC__ 1138 #define __FUNC__ "TSSetOptionsPrefix" 1139 /*@C 1140 TSSetOptionsPrefix - Sets the prefix used for searching for all 1141 TS options in the database. 1142 1143 Collective on TS 1144 1145 Input Parameter: 1146 + ts - The TS context 1147 - prefix - The prefix to prepend to all option names 1148 1149 Notes: 1150 A hyphen (-) must NOT be given at the beginning of the prefix name. 1151 The first character of all runtime options is AUTOMATICALLY the 1152 hyphen. 1153 1154 Contributed by: Matthew Knepley 1155 1156 Level: advanced 1157 1158 .keywords: TS, set, options, prefix, database 1159 1160 .seealso: TSSetFromOptions() 1161 1162 @*/ 1163 int TSSetOptionsPrefix(TS ts, char *prefix) 1164 { 1165 int ierr; 1166 1167 PetscFunctionBegin; 1168 PetscValidHeaderSpecific(ts, TS_COOKIE); 1169 ierr = PetscObjectSetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr); 1170 switch(ts->problem_type) { 1171 case TS_NONLINEAR: 1172 ierr = SNESSetOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr); 1173 break; 1174 case TS_LINEAR: 1175 ierr = SLESSetOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr); 1176 break; 1177 } 1178 PetscFunctionReturn(0); 1179 } 1180 1181 1182 #undef __FUNC__ 1183 #define __FUNC__ "TSAppendOptionsPrefix" 1184 /*@C 1185 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1186 TS options in the database. 1187 1188 Collective on TS 1189 1190 Input Parameter: 1191 + ts - The TS context 1192 - prefix - The prefix to prepend to all option names 1193 1194 Notes: 1195 A hyphen (-) must NOT be given at the beginning of the prefix name. 1196 The first character of all runtime options is AUTOMATICALLY the 1197 hyphen. 1198 1199 Contributed by: Matthew Knepley 1200 1201 Level: advanced 1202 1203 .keywords: TS, append, options, prefix, database 1204 1205 .seealso: TSGetOptionsPrefix() 1206 1207 @*/ 1208 int TSAppendOptionsPrefix(TS ts, char *prefix) 1209 { 1210 int ierr; 1211 1212 PetscFunctionBegin; 1213 PetscValidHeaderSpecific(ts, TS_COOKIE); 1214 ierr = PetscObjectAppendOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr); 1215 switch(ts->problem_type) { 1216 case TS_NONLINEAR: 1217 ierr = SNESAppendOptionsPrefix(ts->snes, prefix);CHKERRQ(ierr); 1218 break; 1219 case TS_LINEAR: 1220 ierr = SLESAppendOptionsPrefix(ts->sles, prefix);CHKERRQ(ierr); 1221 break; 1222 } 1223 PetscFunctionReturn(0); 1224 } 1225 1226 #undef __FUNC__ 1227 #define __FUNC__ "TSGetOptionsPrefix" 1228 /*@C 1229 TSGetOptionsPrefix - Sets the prefix used for searching for all 1230 TS options in the database. 1231 1232 Not Collective 1233 1234 Input Parameter: 1235 . ts - The TS context 1236 1237 Output Parameter: 1238 . prefix - A pointer to the prefix string used 1239 1240 Contributed by: Matthew Knepley 1241 1242 Notes: On the fortran side, the user should pass in a string 'prifix' of 1243 sufficient length to hold the prefix. 1244 1245 Level: intermediate 1246 1247 .keywords: TS, get, options, prefix, database 1248 1249 .seealso: TSAppendOptionsPrefix() 1250 @*/ 1251 int TSGetOptionsPrefix(TS ts, char **prefix) 1252 { 1253 int ierr; 1254 1255 PetscFunctionBegin; 1256 PetscValidHeaderSpecific(ts, TS_COOKIE); 1257 ierr = PetscObjectGetOptionsPrefix((PetscObject) ts, prefix);CHKERRQ(ierr); 1258 PetscFunctionReturn(0); 1259 } 1260 1261 #undef __FUNC__ 1262 #define __FUNC__ "TSGetRHSMatrix" 1263 /*@C 1264 TSGetRHSMatrix - Returns the matrix A at the present timestep. 1265 1266 Not Collective, but parallel objects are returned if TS is parallel 1267 1268 Input Parameter: 1269 . ts - The TS context obtained from TSCreate() 1270 1271 Output Parameters: 1272 + A - The matrix A, where U_t = A(t) U 1273 . M - The preconditioner matrix, usually the same as A 1274 - ctx - User-defined context for matrix evaluation routine 1275 1276 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1277 1278 Contributed by: Matthew Knepley 1279 1280 Level: intermediate 1281 1282 .seealso: TSGetTimeStep(), TSGetTime(), TSGetTimeStepNumber(), TSGetRHSJacobian() 1283 1284 .keywords: TS, timestep, get, matrix 1285 1286 @*/ 1287 int TSGetRHSMatrix(TS ts, Mat *A, Mat *M, void **ctx) 1288 { 1289 PetscFunctionBegin; 1290 PetscValidHeaderSpecific(ts, TS_COOKIE); 1291 if (A) *A = ts->A; 1292 if (M) *M = ts->B; 1293 if (ctx) *ctx = ts->jacP; 1294 PetscFunctionReturn(0); 1295 } 1296 1297 #undef __FUNC__ 1298 #define __FUNC__ "TSGetRHSJacobian" 1299 /*@C 1300 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1301 1302 Not Collective, but parallel objects are returned if TS is parallel 1303 1304 Input Parameter: 1305 . ts - The TS context obtained from TSCreate() 1306 1307 Output Parameters: 1308 + J - The Jacobian J of F, where U_t = F(U,t) 1309 . M - The preconditioner matrix, usually the same as J 1310 - ctx - User-defined context for Jacobian evaluation routine 1311 1312 Notes: You can pass in PETSC_NULL for any return argument you do not need. 1313 1314 Contributed by: Matthew Knepley 1315 1316 Level: intermediate 1317 1318 .seealso: TSGetTimeStep(), TSGetRHSMatrix(), TSGetTime(), TSGetTimeStepNumber() 1319 1320 .keywords: TS, timestep, get, matrix, Jacobian 1321 @*/ 1322 int TSGetRHSJacobian(TS ts, Mat *J, Mat *M, void **ctx) 1323 { 1324 int ierr; 1325 1326 PetscFunctionBegin; 1327 ierr = TSGetRHSMatrix(ts, J, M, ctx);CHKERRQ(ierr); 1328 PetscFunctionReturn(0); 1329 } 1330 1331 /*MC 1332 TSRegisterDynamic - Adds a method to the timestepping solver package. 1333 1334 Synopsis: 1335 1336 TSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(TS)) 1337 1338 Not collective 1339 1340 Input Parameters: 1341 + name_solver - name of a new user-defined solver 1342 . path - path (either absolute or relative) the library containing this solver 1343 . name_create - name of routine to create method context 1344 - routine_create - routine to create method context 1345 1346 Notes: 1347 TSRegisterDynamic() may be called multiple times to add several user-defined solvers. 1348 1349 If dynamic libraries are used, then the fourth input argument (routine_create) 1350 is ignored. 1351 1352 Sample usage: 1353 .vb 1354 TSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 1355 "MySolverCreate",MySolverCreate); 1356 .ve 1357 1358 Then, your solver can be chosen with the procedural interface via 1359 $ TSSetType(ts,"my_solver") 1360 or at runtime via the option 1361 $ -ts_type my_solver 1362 1363 Level: advanced 1364 1365 $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values. 1366 1367 .keywords: TS, register 1368 1369 .seealso: TSRegisterAll(), TSRegisterDestroy() 1370 M*/ 1371 1372 #undef __FUNC__ 1373 #define __FUNC__ "TSRegister" 1374 int TSRegister(char *sname,char *path,char *name,int (*function)(TS)) 1375 { 1376 char fullname[256]; 1377 int ierr; 1378 1379 PetscFunctionBegin; 1380 ierr = FListConcat(path,name,fullname); CHKERRQ(ierr); 1381 ierr = FListAdd(&TSList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr); 1382 PetscFunctionReturn(0); 1383 } 1384