1 2 #include <private/tsimpl.h> /*I "petscts.h" I*/ 3 4 /* Logging support */ 5 PetscClassId TS_CLASSID; 6 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "TSSetTypeFromOptions" 10 /* 11 TSSetTypeFromOptions - Sets the type of ts from user options. 12 13 Collective on TS 14 15 Input Parameter: 16 . ts - The ts 17 18 Level: intermediate 19 20 .keywords: TS, set, options, database, type 21 .seealso: TSSetFromOptions(), TSSetType() 22 */ 23 static PetscErrorCode TSSetTypeFromOptions(TS ts) 24 { 25 PetscBool opt; 26 const char *defaultType; 27 char typeName[256]; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 if (((PetscObject)ts)->type_name) { 32 defaultType = ((PetscObject)ts)->type_name; 33 } else { 34 defaultType = TSEULER; 35 } 36 37 if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 38 ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 39 if (opt) { 40 ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 41 } else { 42 ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 43 } 44 PetscFunctionReturn(0); 45 } 46 47 #undef __FUNCT__ 48 #define __FUNCT__ "TSSetFromOptions" 49 /*@ 50 TSSetFromOptions - Sets various TS parameters from user options. 51 52 Collective on TS 53 54 Input Parameter: 55 . ts - the TS context obtained from TSCreate() 56 57 Options Database Keys: 58 + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 59 . -ts_max_steps maxsteps - maximum number of time-steps to take 60 . -ts_max_time time - maximum time to compute to 61 . -ts_dt dt - initial time step 62 . -ts_monitor - print information at each timestep 63 - -ts_monitor_draw - plot information at each timestep 64 65 Level: beginner 66 67 .keywords: TS, timestep, set, options, database 68 69 .seealso: TSGetType() 70 @*/ 71 PetscErrorCode TSSetFromOptions(TS ts) 72 { 73 PetscReal dt; 74 PetscBool opt,flg; 75 PetscErrorCode ierr; 76 PetscViewer monviewer; 77 char monfilename[PETSC_MAX_PATH_LEN]; 78 SNES snes; 79 80 PetscFunctionBegin; 81 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 82 ierr = PetscOptionsBegin(((PetscObject)ts)->comm, ((PetscObject)ts)->prefix, "Time step options", "TS");CHKERRQ(ierr); 83 84 /* Handle generic TS options */ 85 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 86 ierr = PetscOptionsReal("-ts_max_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 87 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetInitialTime", ts->ptime, &ts->ptime, PETSC_NULL);CHKERRQ(ierr); 88 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetInitialTimeStep",ts->initial_time_step,&dt,&opt);CHKERRQ(ierr); 89 if (opt) { 90 ts->initial_time_step = ts->time_step = dt; 91 } 92 93 /* Monitor options */ 94 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 95 if (flg) { 96 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 97 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 98 } 99 opt = PETSC_FALSE; 100 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 101 if (opt) { 102 ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 103 } 104 opt = PETSC_FALSE; 105 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 106 if (opt) { 107 ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 108 } 109 110 /* Handle TS type options */ 111 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 112 113 /* Handle specific TS options */ 114 if (ts->ops->setfromoptions) { 115 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 116 } 117 118 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 119 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 120 ierr = PetscOptionsEnd();CHKERRQ(ierr); 121 122 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 123 /* Handle subobject options */ 124 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 125 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "TSViewFromOptions" 131 /*@ 132 TSViewFromOptions - This function visualizes the ts based upon user options. 133 134 Collective on TS 135 136 Input Parameter: 137 . ts - The ts 138 139 Level: intermediate 140 141 .keywords: TS, view, options, database 142 .seealso: TSSetFromOptions(), TSView() 143 @*/ 144 PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 145 { 146 PetscViewer viewer; 147 PetscDraw draw; 148 PetscBool opt = PETSC_FALSE; 149 char fileName[PETSC_MAX_PATH_LEN]; 150 PetscErrorCode ierr; 151 152 PetscFunctionBegin; 153 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 154 if (opt && !PetscPreLoadingOn) { 155 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 156 ierr = TSView(ts, viewer);CHKERRQ(ierr); 157 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 158 } 159 opt = PETSC_FALSE; 160 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 161 if (opt) { 162 ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 163 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 164 if (title) { 165 ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 166 } else { 167 ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 168 ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 169 } 170 ierr = TSView(ts, viewer);CHKERRQ(ierr); 171 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 172 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 173 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 174 } 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "TSComputeRHSJacobian" 180 /*@ 181 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 182 set with TSSetRHSJacobian(). 183 184 Collective on TS and Vec 185 186 Input Parameters: 187 + ts - the TS context 188 . t - current timestep 189 - x - input vector 190 191 Output Parameters: 192 + A - Jacobian matrix 193 . B - optional preconditioning matrix 194 - flag - flag indicating matrix structure 195 196 Notes: 197 Most users should not need to explicitly call this routine, as it 198 is used internally within the nonlinear solvers. 199 200 See KSPSetOperators() for important information about setting the 201 flag parameter. 202 203 Level: developer 204 205 .keywords: SNES, compute, Jacobian, matrix 206 207 .seealso: TSSetRHSJacobian(), KSPSetOperators() 208 @*/ 209 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 210 { 211 PetscErrorCode ierr; 212 213 PetscFunctionBegin; 214 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 215 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 216 PetscCheckSameComm(ts,1,X,3); 217 if (ts->ops->rhsjacobian) { 218 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 219 *flg = DIFFERENT_NONZERO_PATTERN; 220 PetscStackPush("TS user Jacobian function"); 221 ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 222 PetscStackPop; 223 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 224 /* make sure user returned a correct Jacobian and preconditioner */ 225 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 226 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 227 } else { 228 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 229 if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 230 *flg = SAME_NONZERO_PATTERN; 231 } 232 PetscFunctionReturn(0); 233 } 234 235 #undef __FUNCT__ 236 #define __FUNCT__ "TSComputeRHSFunction" 237 /*@ 238 TSComputeRHSFunction - Evaluates the right-hand-side function. 239 240 Collective on TS and Vec 241 242 Input Parameters: 243 + ts - the TS context 244 . t - current time 245 - x - state vector 246 247 Output Parameter: 248 . y - right hand side 249 250 Note: 251 Most users should not need to explicitly call this routine, as it 252 is used internally within the nonlinear solvers. 253 254 Level: developer 255 256 .keywords: TS, compute 257 258 .seealso: TSSetRHSFunction(), TSComputeIFunction() 259 @*/ 260 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 261 { 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 266 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 267 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 268 269 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 270 if (ts->ops->rhsfunction) { 271 PetscStackPush("TS user right-hand-side function"); 272 ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 273 PetscStackPop; 274 } else { 275 ierr = VecZeroEntries(y);CHKERRQ(ierr); 276 } 277 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 278 PetscFunctionReturn(0); 279 } 280 281 #undef __FUNCT__ 282 #define __FUNCT__ "TSGetRHSVec_Private" 283 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 284 { 285 PetscErrorCode ierr; 286 287 PetscFunctionBegin; 288 if (!ts->Frhs) { 289 ierr = VecDuplicate(ts->vec_sol,&ts->Frhs);CHKERRQ(ierr); 290 } 291 *Frhs = ts->Frhs; 292 PetscFunctionReturn(0); 293 } 294 295 #undef __FUNCT__ 296 #define __FUNCT__ "TSGetRHSMats_Private" 297 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 298 { 299 PetscErrorCode ierr; 300 Mat A,B; 301 302 PetscFunctionBegin; 303 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 304 if (Arhs) { 305 if (!ts->Arhs) { 306 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 307 } 308 *Arhs = ts->Arhs; 309 } 310 if (Brhs) { 311 if (!ts->Brhs) { 312 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 313 } 314 *Brhs = ts->Brhs; 315 } 316 PetscFunctionReturn(0); 317 } 318 319 #undef __FUNCT__ 320 #define __FUNCT__ "TSComputeIFunction" 321 /*@ 322 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 323 324 Collective on TS and Vec 325 326 Input Parameters: 327 + ts - the TS context 328 . t - current time 329 . X - state vector 330 . Xdot - time derivative of state vector 331 - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 332 333 Output Parameter: 334 . Y - right hand side 335 336 Note: 337 Most users should not need to explicitly call this routine, as it 338 is used internally within the nonlinear solvers. 339 340 If the user did did not write their equations in implicit form, this 341 function recasts them in implicit form. 342 343 Level: developer 344 345 .keywords: TS, compute 346 347 .seealso: TSSetIFunction(), TSComputeRHSFunction() 348 @*/ 349 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 350 { 351 PetscErrorCode ierr; 352 353 PetscFunctionBegin; 354 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 355 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 356 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 357 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 358 359 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 360 if (ts->ops->ifunction) { 361 PetscStackPush("TS user implicit function"); 362 ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 363 PetscStackPop; 364 } 365 if (imex) { 366 if (!ts->ops->ifunction) {ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);} 367 } else { 368 if (!ts->ops->ifunction) { 369 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 370 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 371 } else { 372 Vec Frhs; 373 ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 374 ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 375 ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 376 } 377 } 378 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 379 PetscFunctionReturn(0); 380 } 381 382 #undef __FUNCT__ 383 #define __FUNCT__ "TSComputeIJacobian" 384 /*@ 385 TSComputeIJacobian - Evaluates the Jacobian of the DAE 386 387 Collective on TS and Vec 388 389 Input 390 Input Parameters: 391 + ts - the TS context 392 . t - current timestep 393 . X - state vector 394 . Xdot - time derivative of state vector 395 . shift - shift to apply, see note below 396 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 397 398 Output Parameters: 399 + A - Jacobian matrix 400 . B - optional preconditioning matrix 401 - flag - flag indicating matrix structure 402 403 Notes: 404 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 405 406 dF/dX + shift*dF/dXdot 407 408 Most users should not need to explicitly call this routine, as it 409 is used internally within the nonlinear solvers. 410 411 Level: developer 412 413 .keywords: TS, compute, Jacobian, matrix 414 415 .seealso: TSSetIJacobian() 416 @*/ 417 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 418 { 419 PetscErrorCode ierr; 420 421 PetscFunctionBegin; 422 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 423 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 424 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 425 PetscValidPointer(A,6); 426 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 427 PetscValidPointer(B,7); 428 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 429 PetscValidPointer(flg,8); 430 431 *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 432 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 433 if (ts->ops->ijacobian) { 434 PetscStackPush("TS user implicit Jacobian"); 435 ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 436 PetscStackPop; 437 /* make sure user returned a correct Jacobian and preconditioner */ 438 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 439 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 440 } 441 if (imex) { 442 if (!ts->ops->ijacobian) { /* system was written as Xdot = F(t,X) */ 443 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 444 ierr = MatShift(*A,1.0);CHKERRQ(ierr); 445 if (*A != *B) { 446 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 447 ierr = MatShift(*B,shift);CHKERRQ(ierr); 448 } 449 *flg = SAME_PRECONDITIONER; 450 } 451 } else { 452 if (!ts->ops->ijacobian) { 453 ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 454 ierr = MatScale(*A,-1);CHKERRQ(ierr); 455 ierr = MatShift(*A,shift);CHKERRQ(ierr); 456 if (*A != *B) { 457 ierr = MatScale(*B,-1);CHKERRQ(ierr); 458 ierr = MatShift(*B,shift);CHKERRQ(ierr); 459 } 460 } else if (ts->ops->rhsjacobian) { 461 Mat Arhs,Brhs; 462 MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 463 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 464 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 465 axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 466 ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 467 if (*A != *B) { 468 ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 469 } 470 *flg = PetscMin(*flg,flg2); 471 } 472 } 473 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 474 PetscFunctionReturn(0); 475 } 476 477 #undef __FUNCT__ 478 #define __FUNCT__ "TSSetRHSFunction" 479 /*@C 480 TSSetRHSFunction - Sets the routine for evaluating the function, 481 F(t,u), where U_t = F(t,u). 482 483 Logically Collective on TS 484 485 Input Parameters: 486 + ts - the TS context obtained from TSCreate() 487 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 488 . f - routine for evaluating the right-hand-side function 489 - ctx - [optional] user-defined context for private data for the 490 function evaluation routine (may be PETSC_NULL) 491 492 Calling sequence of func: 493 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 494 495 + t - current timestep 496 . u - input vector 497 . F - function vector 498 - ctx - [optional] user-defined function context 499 500 Important: 501 The user MUST call either this routine or TSSetMatrices(). 502 503 Level: beginner 504 505 .keywords: TS, timestep, set, right-hand-side, function 506 507 .seealso: TSSetMatrices() 508 @*/ 509 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 510 { 511 PetscErrorCode ierr; 512 SNES snes; 513 514 PetscFunctionBegin; 515 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 516 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 517 if (f) ts->ops->rhsfunction = f; 518 if (ctx) ts->funP = ctx; 519 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 520 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 521 PetscFunctionReturn(0); 522 } 523 524 #undef __FUNCT__ 525 #define __FUNCT__ "TSSetRHSJacobian" 526 /*@C 527 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 528 where U_t = F(U,t), as well as the location to store the matrix. 529 Use TSSetMatrices() for linear problems. 530 531 Logically Collective on TS 532 533 Input Parameters: 534 + ts - the TS context obtained from TSCreate() 535 . A - Jacobian matrix 536 . B - preconditioner matrix (usually same as A) 537 . f - the Jacobian evaluation routine 538 - ctx - [optional] user-defined context for private data for the 539 Jacobian evaluation routine (may be PETSC_NULL) 540 541 Calling sequence of func: 542 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 543 544 + t - current timestep 545 . u - input vector 546 . A - matrix A, where U_t = A(t)u 547 . B - preconditioner matrix, usually the same as A 548 . flag - flag indicating information about the preconditioner matrix 549 structure (same as flag in KSPSetOperators()) 550 - ctx - [optional] user-defined context for matrix evaluation routine 551 552 Notes: 553 See KSPSetOperators() for important information about setting the flag 554 output parameter in the routine func(). Be sure to read this information! 555 556 The routine func() takes Mat * as the matrix arguments rather than Mat. 557 This allows the matrix evaluation routine to replace A and/or B with a 558 completely new matrix structure (not just different matrix elements) 559 when appropriate, for instance, if the nonzero structure is changing 560 throughout the global iterations. 561 562 Level: beginner 563 564 .keywords: TS, timestep, set, right-hand-side, Jacobian 565 566 .seealso: TSDefaultComputeJacobianColor(), 567 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 568 569 @*/ 570 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 571 { 572 PetscErrorCode ierr; 573 SNES snes; 574 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 577 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 578 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 579 if (A) PetscCheckSameComm(ts,1,A,2); 580 if (B) PetscCheckSameComm(ts,1,B,3); 581 582 if (f) ts->ops->rhsjacobian = f; 583 if (ctx) ts->jacP = ctx; 584 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 585 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 586 PetscFunctionReturn(0); 587 } 588 589 590 #undef __FUNCT__ 591 #define __FUNCT__ "TSSetIFunction" 592 /*@C 593 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 594 595 Logically Collective on TS 596 597 Input Parameters: 598 + ts - the TS context obtained from TSCreate() 599 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 600 . f - the function evaluation routine 601 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 602 603 Calling sequence of f: 604 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 605 606 + t - time at step/stage being solved 607 . u - state vector 608 . u_t - time derivative of state vector 609 . F - function vector 610 - ctx - [optional] user-defined context for matrix evaluation routine 611 612 Important: 613 The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 614 615 Level: beginner 616 617 .keywords: TS, timestep, set, DAE, Jacobian 618 619 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 620 @*/ 621 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 622 { 623 PetscErrorCode ierr; 624 SNES snes; 625 626 PetscFunctionBegin; 627 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 628 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 629 if (f) ts->ops->ifunction = f; 630 if (ctx) ts->funP = ctx; 631 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 632 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 633 PetscFunctionReturn(0); 634 } 635 636 #undef __FUNCT__ 637 #define __FUNCT__ "TSGetIFunction" 638 /*@C 639 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 640 641 Not Collective 642 643 Input Parameter: 644 . ts - the TS context 645 646 Output Parameter: 647 + r - vector to hold residual (or PETSC_NULL) 648 . func - the function to compute residual (or PETSC_NULL) 649 - ctx - the function context (or PETSC_NULL) 650 651 Level: advanced 652 653 .keywords: TS, nonlinear, get, function 654 655 .seealso: TSSetIFunction(), SNESGetFunction() 656 @*/ 657 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 658 { 659 PetscErrorCode ierr; 660 SNES snes; 661 662 PetscFunctionBegin; 663 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 664 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 665 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 666 if (func) *func = ts->ops->ifunction; 667 if (ctx) *ctx = ts->funP; 668 PetscFunctionReturn(0); 669 } 670 671 #undef __FUNCT__ 672 #define __FUNCT__ "TSGetRHSFunction" 673 /*@C 674 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 675 676 Not Collective 677 678 Input Parameter: 679 . ts - the TS context 680 681 Output Parameter: 682 + r - vector to hold computed right hand side (or PETSC_NULL) 683 . func - the function to compute right hand side (or PETSC_NULL) 684 - ctx - the function context (or PETSC_NULL) 685 686 Level: advanced 687 688 .keywords: TS, nonlinear, get, function 689 690 .seealso: TSSetRhsfunction(), SNESGetFunction() 691 @*/ 692 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 693 { 694 PetscErrorCode ierr; 695 SNES snes; 696 697 PetscFunctionBegin; 698 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 699 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 700 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 701 if (func) *func = ts->ops->rhsfunction; 702 if (ctx) *ctx = ts->funP; 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "TSSetIJacobian" 708 /*@C 709 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 710 you provided with TSSetIFunction(). 711 712 Logically Collective on TS 713 714 Input Parameters: 715 + ts - the TS context obtained from TSCreate() 716 . A - Jacobian matrix 717 . B - preconditioning matrix for A (may be same as A) 718 . f - the Jacobian evaluation routine 719 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 720 721 Calling sequence of f: 722 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 723 724 + t - time at step/stage being solved 725 . U - state vector 726 . U_t - time derivative of state vector 727 . a - shift 728 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 729 . B - preconditioning matrix for A, may be same as A 730 . flag - flag indicating information about the preconditioner matrix 731 structure (same as flag in KSPSetOperators()) 732 - ctx - [optional] user-defined context for matrix evaluation routine 733 734 Notes: 735 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 736 737 The matrix dF/dU + a*dF/dU_t you provide turns out to be 738 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 739 The time integrator internally approximates U_t by W+a*U where the positive "shift" 740 a and vector W depend on the integration method, step size, and past states. For example with 741 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 742 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 743 744 Level: beginner 745 746 .keywords: TS, timestep, DAE, Jacobian 747 748 .seealso: TSSetIFunction(), TSSetRHSJacobian() 749 750 @*/ 751 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 752 { 753 PetscErrorCode ierr; 754 SNES snes; 755 756 PetscFunctionBegin; 757 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 758 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 759 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 760 if (A) PetscCheckSameComm(ts,1,A,2); 761 if (B) PetscCheckSameComm(ts,1,B,3); 762 if (f) ts->ops->ijacobian = f; 763 if (ctx) ts->jacP = ctx; 764 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 765 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 766 PetscFunctionReturn(0); 767 } 768 769 #undef __FUNCT__ 770 #define __FUNCT__ "TSView" 771 /*@C 772 TSView - Prints the TS data structure. 773 774 Collective on TS 775 776 Input Parameters: 777 + ts - the TS context obtained from TSCreate() 778 - viewer - visualization context 779 780 Options Database Key: 781 . -ts_view - calls TSView() at end of TSStep() 782 783 Notes: 784 The available visualization contexts include 785 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 786 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 787 output where only the first processor opens 788 the file. All other processors send their 789 data to the first processor to print. 790 791 The user can open an alternative visualization context with 792 PetscViewerASCIIOpen() - output to a specified file. 793 794 Level: beginner 795 796 .keywords: TS, timestep, view 797 798 .seealso: PetscViewerASCIIOpen() 799 @*/ 800 PetscErrorCode TSView(TS ts,PetscViewer viewer) 801 { 802 PetscErrorCode ierr; 803 const TSType type; 804 PetscBool iascii,isstring,isundials; 805 806 PetscFunctionBegin; 807 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 808 if (!viewer) { 809 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 810 } 811 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 812 PetscCheckSameComm(ts,1,viewer,2); 813 814 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 815 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 816 if (iascii) { 817 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 818 if (ts->ops->view) { 819 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 820 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 821 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 822 } 823 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 824 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 825 if (ts->problem_type == TS_NONLINEAR) { 826 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 827 } 828 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 829 } else if (isstring) { 830 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 831 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 832 } 833 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 834 ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 835 if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 836 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 837 PetscFunctionReturn(0); 838 } 839 840 841 #undef __FUNCT__ 842 #define __FUNCT__ "TSSetApplicationContext" 843 /*@ 844 TSSetApplicationContext - Sets an optional user-defined context for 845 the timesteppers. 846 847 Logically Collective on TS 848 849 Input Parameters: 850 + ts - the TS context obtained from TSCreate() 851 - usrP - optional user context 852 853 Level: intermediate 854 855 .keywords: TS, timestep, set, application, context 856 857 .seealso: TSGetApplicationContext() 858 @*/ 859 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 860 { 861 PetscFunctionBegin; 862 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 863 ts->user = usrP; 864 PetscFunctionReturn(0); 865 } 866 867 #undef __FUNCT__ 868 #define __FUNCT__ "TSGetApplicationContext" 869 /*@ 870 TSGetApplicationContext - Gets the user-defined context for the 871 timestepper. 872 873 Not Collective 874 875 Input Parameter: 876 . ts - the TS context obtained from TSCreate() 877 878 Output Parameter: 879 . usrP - user context 880 881 Level: intermediate 882 883 .keywords: TS, timestep, get, application, context 884 885 .seealso: TSSetApplicationContext() 886 @*/ 887 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 888 { 889 PetscFunctionBegin; 890 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 891 *(void**)usrP = ts->user; 892 PetscFunctionReturn(0); 893 } 894 895 #undef __FUNCT__ 896 #define __FUNCT__ "TSGetTimeStepNumber" 897 /*@ 898 TSGetTimeStepNumber - Gets the current number of timesteps. 899 900 Not Collective 901 902 Input Parameter: 903 . ts - the TS context obtained from TSCreate() 904 905 Output Parameter: 906 . iter - number steps so far 907 908 Level: intermediate 909 910 .keywords: TS, timestep, get, iteration, number 911 @*/ 912 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 913 { 914 PetscFunctionBegin; 915 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 916 PetscValidIntPointer(iter,2); 917 *iter = ts->steps; 918 PetscFunctionReturn(0); 919 } 920 921 #undef __FUNCT__ 922 #define __FUNCT__ "TSSetInitialTimeStep" 923 /*@ 924 TSSetInitialTimeStep - Sets the initial timestep to be used, 925 as well as the initial time. 926 927 Logically Collective on TS 928 929 Input Parameters: 930 + ts - the TS context obtained from TSCreate() 931 . initial_time - the initial time 932 - time_step - the size of the timestep 933 934 Level: intermediate 935 936 .seealso: TSSetTimeStep(), TSGetTimeStep() 937 938 .keywords: TS, set, initial, timestep 939 @*/ 940 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 941 { 942 PetscFunctionBegin; 943 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 944 ts->time_step = time_step; 945 ts->initial_time_step = time_step; 946 ts->ptime = initial_time; 947 PetscFunctionReturn(0); 948 } 949 950 #undef __FUNCT__ 951 #define __FUNCT__ "TSSetTimeStep" 952 /*@ 953 TSSetTimeStep - Allows one to reset the timestep at any time, 954 useful for simple pseudo-timestepping codes. 955 956 Logically Collective on TS 957 958 Input Parameters: 959 + ts - the TS context obtained from TSCreate() 960 - time_step - the size of the timestep 961 962 Level: intermediate 963 964 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 965 966 .keywords: TS, set, timestep 967 @*/ 968 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 969 { 970 PetscFunctionBegin; 971 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 972 PetscValidLogicalCollectiveReal(ts,time_step,2); 973 ts->time_step = time_step; 974 PetscFunctionReturn(0); 975 } 976 977 #undef __FUNCT__ 978 #define __FUNCT__ "TSGetTimeStep" 979 /*@ 980 TSGetTimeStep - Gets the current timestep size. 981 982 Not Collective 983 984 Input Parameter: 985 . ts - the TS context obtained from TSCreate() 986 987 Output Parameter: 988 . dt - the current timestep size 989 990 Level: intermediate 991 992 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 993 994 .keywords: TS, get, timestep 995 @*/ 996 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 997 { 998 PetscFunctionBegin; 999 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1000 PetscValidDoublePointer(dt,2); 1001 *dt = ts->time_step; 1002 PetscFunctionReturn(0); 1003 } 1004 1005 #undef __FUNCT__ 1006 #define __FUNCT__ "TSGetSolution" 1007 /*@ 1008 TSGetSolution - Returns the solution at the present timestep. It 1009 is valid to call this routine inside the function that you are evaluating 1010 in order to move to the new timestep. This vector not changed until 1011 the solution at the next timestep has been calculated. 1012 1013 Not Collective, but Vec returned is parallel if TS is parallel 1014 1015 Input Parameter: 1016 . ts - the TS context obtained from TSCreate() 1017 1018 Output Parameter: 1019 . v - the vector containing the solution 1020 1021 Level: intermediate 1022 1023 .seealso: TSGetTimeStep() 1024 1025 .keywords: TS, timestep, get, solution 1026 @*/ 1027 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1028 { 1029 PetscFunctionBegin; 1030 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1031 PetscValidPointer(v,2); 1032 *v = ts->vec_sol; 1033 PetscFunctionReturn(0); 1034 } 1035 1036 /* ----- Routines to initialize and destroy a timestepper ---- */ 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "TSSetProblemType" 1039 /*@ 1040 TSSetProblemType - Sets the type of problem to be solved. 1041 1042 Not collective 1043 1044 Input Parameters: 1045 + ts - The TS 1046 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1047 .vb 1048 U_t = A U 1049 U_t = A(t) U 1050 U_t = F(t,U) 1051 .ve 1052 1053 Level: beginner 1054 1055 .keywords: TS, problem type 1056 .seealso: TSSetUp(), TSProblemType, TS 1057 @*/ 1058 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1059 { 1060 PetscErrorCode ierr; 1061 1062 PetscFunctionBegin; 1063 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1064 ts->problem_type = type; 1065 if (type == TS_LINEAR) { 1066 SNES snes; 1067 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1068 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 1069 } 1070 PetscFunctionReturn(0); 1071 } 1072 1073 #undef __FUNCT__ 1074 #define __FUNCT__ "TSGetProblemType" 1075 /*@C 1076 TSGetProblemType - Gets the type of problem to be solved. 1077 1078 Not collective 1079 1080 Input Parameter: 1081 . ts - The TS 1082 1083 Output Parameter: 1084 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1085 .vb 1086 M U_t = A U 1087 M(t) U_t = A(t) U 1088 U_t = F(t,U) 1089 .ve 1090 1091 Level: beginner 1092 1093 .keywords: TS, problem type 1094 .seealso: TSSetUp(), TSProblemType, TS 1095 @*/ 1096 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1097 { 1098 PetscFunctionBegin; 1099 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1100 PetscValidIntPointer(type,2); 1101 *type = ts->problem_type; 1102 PetscFunctionReturn(0); 1103 } 1104 1105 #undef __FUNCT__ 1106 #define __FUNCT__ "TSSetUp" 1107 /*@ 1108 TSSetUp - Sets up the internal data structures for the later use 1109 of a timestepper. 1110 1111 Collective on TS 1112 1113 Input Parameter: 1114 . ts - the TS context obtained from TSCreate() 1115 1116 Notes: 1117 For basic use of the TS solvers the user need not explicitly call 1118 TSSetUp(), since these actions will automatically occur during 1119 the call to TSStep(). However, if one wishes to control this 1120 phase separately, TSSetUp() should be called after TSCreate() 1121 and optional routines of the form TSSetXXX(), but before TSStep(). 1122 1123 Level: advanced 1124 1125 .keywords: TS, timestep, setup 1126 1127 .seealso: TSCreate(), TSStep(), TSDestroy() 1128 @*/ 1129 PetscErrorCode TSSetUp(TS ts) 1130 { 1131 PetscErrorCode ierr; 1132 1133 PetscFunctionBegin; 1134 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1135 if (ts->setupcalled) PetscFunctionReturn(0); 1136 1137 if (!((PetscObject)ts)->type_name) { 1138 ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1139 } 1140 1141 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1142 1143 if (ts->ops->setup) { 1144 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1145 } 1146 1147 ts->setupcalled = PETSC_TRUE; 1148 PetscFunctionReturn(0); 1149 } 1150 1151 #undef __FUNCT__ 1152 #define __FUNCT__ "TSReset" 1153 /*@ 1154 TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1155 1156 Collective on TS 1157 1158 Input Parameter: 1159 . ts - the TS context obtained from TSCreate() 1160 1161 Level: beginner 1162 1163 .keywords: TS, timestep, reset 1164 1165 .seealso: TSCreate(), TSSetup(), TSDestroy() 1166 @*/ 1167 PetscErrorCode TSReset(TS ts) 1168 { 1169 PetscErrorCode ierr; 1170 1171 PetscFunctionBegin; 1172 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1173 if (ts->ops->reset) { 1174 ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1175 } 1176 if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1177 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 1178 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 1179 ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr); 1180 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1181 if (ts->work) {ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr);} 1182 ts->setupcalled = PETSC_FALSE; 1183 PetscFunctionReturn(0); 1184 } 1185 1186 #undef __FUNCT__ 1187 #define __FUNCT__ "TSDestroy" 1188 /*@ 1189 TSDestroy - Destroys the timestepper context that was created 1190 with TSCreate(). 1191 1192 Collective on TS 1193 1194 Input Parameter: 1195 . ts - the TS context obtained from TSCreate() 1196 1197 Level: beginner 1198 1199 .keywords: TS, timestepper, destroy 1200 1201 .seealso: TSCreate(), TSSetUp(), TSSolve() 1202 @*/ 1203 PetscErrorCode TSDestroy(TS *ts) 1204 { 1205 PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 if (!*ts) PetscFunctionReturn(0); 1209 PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 1210 if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1211 1212 ierr = TSReset((*ts));CHKERRQ(ierr); 1213 1214 /* if memory was published with AMS then destroy it */ 1215 ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 1216 if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 1217 1218 ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 1219 ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 1220 ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 1221 1222 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1223 PetscFunctionReturn(0); 1224 } 1225 1226 #undef __FUNCT__ 1227 #define __FUNCT__ "TSGetSNES" 1228 /*@ 1229 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1230 a TS (timestepper) context. Valid only for nonlinear problems. 1231 1232 Not Collective, but SNES is parallel if TS is parallel 1233 1234 Input Parameter: 1235 . ts - the TS context obtained from TSCreate() 1236 1237 Output Parameter: 1238 . snes - the nonlinear solver context 1239 1240 Notes: 1241 The user can then directly manipulate the SNES context to set various 1242 options, etc. Likewise, the user can then extract and manipulate the 1243 KSP, KSP, and PC contexts as well. 1244 1245 TSGetSNES() does not work for integrators that do not use SNES; in 1246 this case TSGetSNES() returns PETSC_NULL in snes. 1247 1248 Level: beginner 1249 1250 .keywords: timestep, get, SNES 1251 @*/ 1252 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1253 { 1254 PetscErrorCode ierr; 1255 1256 PetscFunctionBegin; 1257 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1258 PetscValidPointer(snes,2); 1259 if (!ts->snes) { 1260 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1261 ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1262 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1263 if (ts->problem_type == TS_LINEAR) { 1264 ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr); 1265 } 1266 } 1267 *snes = ts->snes; 1268 PetscFunctionReturn(0); 1269 } 1270 1271 #undef __FUNCT__ 1272 #define __FUNCT__ "TSGetKSP" 1273 /*@ 1274 TSGetKSP - Returns the KSP (linear solver) associated with 1275 a TS (timestepper) context. 1276 1277 Not Collective, but KSP is parallel if TS is parallel 1278 1279 Input Parameter: 1280 . ts - the TS context obtained from TSCreate() 1281 1282 Output Parameter: 1283 . ksp - the nonlinear solver context 1284 1285 Notes: 1286 The user can then directly manipulate the KSP context to set various 1287 options, etc. Likewise, the user can then extract and manipulate the 1288 KSP and PC contexts as well. 1289 1290 TSGetKSP() does not work for integrators that do not use KSP; 1291 in this case TSGetKSP() returns PETSC_NULL in ksp. 1292 1293 Level: beginner 1294 1295 .keywords: timestep, get, KSP 1296 @*/ 1297 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1298 { 1299 PetscErrorCode ierr; 1300 SNES snes; 1301 1302 PetscFunctionBegin; 1303 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1304 PetscValidPointer(ksp,2); 1305 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1306 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1307 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1308 ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1309 PetscFunctionReturn(0); 1310 } 1311 1312 /* ----------- Routines to set solver parameters ---------- */ 1313 1314 #undef __FUNCT__ 1315 #define __FUNCT__ "TSGetDuration" 1316 /*@ 1317 TSGetDuration - Gets the maximum number of timesteps to use and 1318 maximum time for iteration. 1319 1320 Not Collective 1321 1322 Input Parameters: 1323 + ts - the TS context obtained from TSCreate() 1324 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1325 - maxtime - final time to iterate to, or PETSC_NULL 1326 1327 Level: intermediate 1328 1329 .keywords: TS, timestep, get, maximum, iterations, time 1330 @*/ 1331 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1332 { 1333 PetscFunctionBegin; 1334 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1335 if (maxsteps) { 1336 PetscValidIntPointer(maxsteps,2); 1337 *maxsteps = ts->max_steps; 1338 } 1339 if (maxtime ) { 1340 PetscValidScalarPointer(maxtime,3); 1341 *maxtime = ts->max_time; 1342 } 1343 PetscFunctionReturn(0); 1344 } 1345 1346 #undef __FUNCT__ 1347 #define __FUNCT__ "TSSetDuration" 1348 /*@ 1349 TSSetDuration - Sets the maximum number of timesteps to use and 1350 maximum time for iteration. 1351 1352 Logically Collective on TS 1353 1354 Input Parameters: 1355 + ts - the TS context obtained from TSCreate() 1356 . maxsteps - maximum number of iterations to use 1357 - maxtime - final time to iterate to 1358 1359 Options Database Keys: 1360 . -ts_max_steps <maxsteps> - Sets maxsteps 1361 . -ts_max_time <maxtime> - Sets maxtime 1362 1363 Notes: 1364 The default maximum number of iterations is 5000. Default time is 5.0 1365 1366 Level: intermediate 1367 1368 .keywords: TS, timestep, set, maximum, iterations 1369 @*/ 1370 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1371 { 1372 PetscFunctionBegin; 1373 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1374 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1375 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1376 ts->max_steps = maxsteps; 1377 ts->max_time = maxtime; 1378 PetscFunctionReturn(0); 1379 } 1380 1381 #undef __FUNCT__ 1382 #define __FUNCT__ "TSSetSolution" 1383 /*@ 1384 TSSetSolution - Sets the initial solution vector 1385 for use by the TS routines. 1386 1387 Logically Collective on TS and Vec 1388 1389 Input Parameters: 1390 + ts - the TS context obtained from TSCreate() 1391 - x - the solution vector 1392 1393 Level: beginner 1394 1395 .keywords: TS, timestep, set, solution, initial conditions 1396 @*/ 1397 PetscErrorCode TSSetSolution(TS ts,Vec x) 1398 { 1399 PetscErrorCode ierr; 1400 1401 PetscFunctionBegin; 1402 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1403 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1404 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1405 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1406 ts->vec_sol = x; 1407 PetscFunctionReturn(0); 1408 } 1409 1410 #undef __FUNCT__ 1411 #define __FUNCT__ "TSSetPreStep" 1412 /*@C 1413 TSSetPreStep - Sets the general-purpose function 1414 called once at the beginning of each time step. 1415 1416 Logically Collective on TS 1417 1418 Input Parameters: 1419 + ts - The TS context obtained from TSCreate() 1420 - func - The function 1421 1422 Calling sequence of func: 1423 . func (TS ts); 1424 1425 Level: intermediate 1426 1427 .keywords: TS, timestep 1428 @*/ 1429 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1430 { 1431 PetscFunctionBegin; 1432 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1433 ts->ops->prestep = func; 1434 PetscFunctionReturn(0); 1435 } 1436 1437 #undef __FUNCT__ 1438 #define __FUNCT__ "TSPreStep" 1439 /*@C 1440 TSPreStep - Runs the user-defined pre-step function. 1441 1442 Collective on TS 1443 1444 Input Parameters: 1445 . ts - The TS context obtained from TSCreate() 1446 1447 Notes: 1448 TSPreStep() is typically used within time stepping implementations, 1449 so most users would not generally call this routine themselves. 1450 1451 Level: developer 1452 1453 .keywords: TS, timestep 1454 @*/ 1455 PetscErrorCode TSPreStep(TS ts) 1456 { 1457 PetscErrorCode ierr; 1458 1459 PetscFunctionBegin; 1460 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1461 if (ts->ops->prestep) { 1462 PetscStackPush("TS PreStep function"); 1463 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1464 PetscStackPop; 1465 } 1466 PetscFunctionReturn(0); 1467 } 1468 1469 #undef __FUNCT__ 1470 #define __FUNCT__ "TSSetPostStep" 1471 /*@C 1472 TSSetPostStep - Sets the general-purpose function 1473 called once at the end of each time step. 1474 1475 Logically Collective on TS 1476 1477 Input Parameters: 1478 + ts - The TS context obtained from TSCreate() 1479 - func - The function 1480 1481 Calling sequence of func: 1482 . func (TS ts); 1483 1484 Level: intermediate 1485 1486 .keywords: TS, timestep 1487 @*/ 1488 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1489 { 1490 PetscFunctionBegin; 1491 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1492 ts->ops->poststep = func; 1493 PetscFunctionReturn(0); 1494 } 1495 1496 #undef __FUNCT__ 1497 #define __FUNCT__ "TSPostStep" 1498 /*@C 1499 TSPostStep - Runs the user-defined post-step function. 1500 1501 Collective on TS 1502 1503 Input Parameters: 1504 . ts - The TS context obtained from TSCreate() 1505 1506 Notes: 1507 TSPostStep() is typically used within time stepping implementations, 1508 so most users would not generally call this routine themselves. 1509 1510 Level: developer 1511 1512 .keywords: TS, timestep 1513 @*/ 1514 PetscErrorCode TSPostStep(TS ts) 1515 { 1516 PetscErrorCode ierr; 1517 1518 PetscFunctionBegin; 1519 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1520 if (ts->ops->poststep) { 1521 PetscStackPush("TS PostStep function"); 1522 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1523 PetscStackPop; 1524 } 1525 PetscFunctionReturn(0); 1526 } 1527 1528 /* ------------ Routines to set performance monitoring options ----------- */ 1529 1530 #undef __FUNCT__ 1531 #define __FUNCT__ "TSMonitorSet" 1532 /*@C 1533 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1534 timestep to display the iteration's progress. 1535 1536 Logically Collective on TS 1537 1538 Input Parameters: 1539 + ts - the TS context obtained from TSCreate() 1540 . func - monitoring routine 1541 . mctx - [optional] user-defined context for private data for the 1542 monitor routine (use PETSC_NULL if no context is desired) 1543 - monitordestroy - [optional] routine that frees monitor context 1544 (may be PETSC_NULL) 1545 1546 Calling sequence of func: 1547 $ int func(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1548 1549 + ts - the TS context 1550 . steps - iteration number 1551 . time - current time 1552 . x - current iterate 1553 - mctx - [optional] monitoring context 1554 1555 Notes: 1556 This routine adds an additional monitor to the list of monitors that 1557 already has been loaded. 1558 1559 Fortran notes: Only a single monitor function can be set for each TS object 1560 1561 Level: intermediate 1562 1563 .keywords: TS, timestep, set, monitor 1564 1565 .seealso: TSMonitorDefault(), TSMonitorCancel() 1566 @*/ 1567 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1568 { 1569 PetscFunctionBegin; 1570 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1571 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1572 ts->monitor[ts->numbermonitors] = monitor; 1573 ts->mdestroy[ts->numbermonitors] = mdestroy; 1574 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1575 PetscFunctionReturn(0); 1576 } 1577 1578 #undef __FUNCT__ 1579 #define __FUNCT__ "TSMonitorCancel" 1580 /*@C 1581 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1582 1583 Logically Collective on TS 1584 1585 Input Parameters: 1586 . ts - the TS context obtained from TSCreate() 1587 1588 Notes: 1589 There is no way to remove a single, specific monitor. 1590 1591 Level: intermediate 1592 1593 .keywords: TS, timestep, set, monitor 1594 1595 .seealso: TSMonitorDefault(), TSMonitorSet() 1596 @*/ 1597 PetscErrorCode TSMonitorCancel(TS ts) 1598 { 1599 PetscErrorCode ierr; 1600 PetscInt i; 1601 1602 PetscFunctionBegin; 1603 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1604 for (i=0; i<ts->numbermonitors; i++) { 1605 if (ts->mdestroy[i]) { 1606 ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1607 } 1608 } 1609 ts->numbermonitors = 0; 1610 PetscFunctionReturn(0); 1611 } 1612 1613 #undef __FUNCT__ 1614 #define __FUNCT__ "TSMonitorDefault" 1615 /*@ 1616 TSMonitorDefault - Sets the Default monitor 1617 1618 Level: intermediate 1619 1620 .keywords: TS, set, monitor 1621 1622 .seealso: TSMonitorDefault(), TSMonitorSet() 1623 @*/ 1624 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1625 { 1626 PetscErrorCode ierr; 1627 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1628 1629 PetscFunctionBegin; 1630 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1631 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %G time %G\n",step,ts->time_step,ptime);CHKERRQ(ierr); 1632 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1633 PetscFunctionReturn(0); 1634 } 1635 1636 #undef __FUNCT__ 1637 #define __FUNCT__ "TSStep" 1638 /*@ 1639 TSStep - Steps the requested number of timesteps. 1640 1641 Collective on TS 1642 1643 Input Parameter: 1644 . ts - the TS context obtained from TSCreate() 1645 1646 Output Parameters: 1647 + steps - number of iterations until termination 1648 - ptime - time until termination 1649 1650 Level: beginner 1651 1652 .keywords: TS, timestep, solve 1653 1654 .seealso: TSCreate(), TSSetUp(), TSDestroy() 1655 @*/ 1656 PetscErrorCode TSStep(TS ts,PetscInt *steps,PetscReal *ptime) 1657 { 1658 PetscErrorCode ierr; 1659 1660 PetscFunctionBegin; 1661 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1662 1663 ierr = TSSetUp(ts);CHKERRQ(ierr); 1664 1665 ierr = PetscLogEventBegin(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1666 ierr = (*ts->ops->step)(ts, steps, ptime);CHKERRQ(ierr); 1667 ierr = PetscLogEventEnd(TS_Step, ts, 0, 0, 0);CHKERRQ(ierr); 1668 1669 if (!PetscPreLoadingOn) { 1670 ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1671 } 1672 PetscFunctionReturn(0); 1673 } 1674 1675 #undef __FUNCT__ 1676 #define __FUNCT__ "TSSolve" 1677 /*@ 1678 TSSolve - Steps the requested number of timesteps. 1679 1680 Collective on TS 1681 1682 Input Parameter: 1683 + ts - the TS context obtained from TSCreate() 1684 - x - the solution vector, or PETSC_NULL if it was set with TSSetSolution() 1685 1686 Level: beginner 1687 1688 .keywords: TS, timestep, solve 1689 1690 .seealso: TSCreate(), TSSetSolution(), TSStep() 1691 @*/ 1692 PetscErrorCode TSSolve(TS ts, Vec x) 1693 { 1694 PetscInt steps; 1695 PetscReal ptime; 1696 PetscErrorCode ierr; 1697 1698 PetscFunctionBegin; 1699 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1700 /* set solution vector if provided */ 1701 if (x) { ierr = TSSetSolution(ts, x); CHKERRQ(ierr); } 1702 /* reset time step and iteration counters */ 1703 ts->steps = 0; ts->linear_its = 0; ts->nonlinear_its = 0; 1704 /* steps the requested number of timesteps. */ 1705 ierr = TSStep(ts, &steps, &ptime);CHKERRQ(ierr); 1706 PetscFunctionReturn(0); 1707 } 1708 1709 #undef __FUNCT__ 1710 #define __FUNCT__ "TSMonitor" 1711 /* 1712 Runs the user provided monitor routines, if they exists. 1713 */ 1714 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1715 { 1716 PetscErrorCode ierr; 1717 PetscInt i,n = ts->numbermonitors; 1718 1719 PetscFunctionBegin; 1720 for (i=0; i<n; i++) { 1721 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1722 } 1723 PetscFunctionReturn(0); 1724 } 1725 1726 /* ------------------------------------------------------------------------*/ 1727 1728 #undef __FUNCT__ 1729 #define __FUNCT__ "TSMonitorLGCreate" 1730 /*@C 1731 TSMonitorLGCreate - Creates a line graph context for use with 1732 TS to monitor convergence of preconditioned residual norms. 1733 1734 Collective on TS 1735 1736 Input Parameters: 1737 + host - the X display to open, or null for the local machine 1738 . label - the title to put in the title bar 1739 . x, y - the screen coordinates of the upper left coordinate of the window 1740 - m, n - the screen width and height in pixels 1741 1742 Output Parameter: 1743 . draw - the drawing context 1744 1745 Options Database Key: 1746 . -ts_monitor_draw - automatically sets line graph monitor 1747 1748 Notes: 1749 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1750 1751 Level: intermediate 1752 1753 .keywords: TS, monitor, line graph, residual, seealso 1754 1755 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1756 1757 @*/ 1758 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1759 { 1760 PetscDraw win; 1761 PetscErrorCode ierr; 1762 1763 PetscFunctionBegin; 1764 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1765 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1766 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1767 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1768 1769 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1770 PetscFunctionReturn(0); 1771 } 1772 1773 #undef __FUNCT__ 1774 #define __FUNCT__ "TSMonitorLG" 1775 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1776 { 1777 PetscDrawLG lg = (PetscDrawLG) monctx; 1778 PetscReal x,y = ptime; 1779 PetscErrorCode ierr; 1780 1781 PetscFunctionBegin; 1782 if (!monctx) { 1783 MPI_Comm comm; 1784 PetscViewer viewer; 1785 1786 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1787 viewer = PETSC_VIEWER_DRAW_(comm); 1788 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1789 } 1790 1791 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1792 x = (PetscReal)n; 1793 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1794 if (n < 20 || (n % 5)) { 1795 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1796 } 1797 PetscFunctionReturn(0); 1798 } 1799 1800 #undef __FUNCT__ 1801 #define __FUNCT__ "TSMonitorLGDestroy" 1802 /*@C 1803 TSMonitorLGDestroy - Destroys a line graph context that was created 1804 with TSMonitorLGCreate(). 1805 1806 Collective on PetscDrawLG 1807 1808 Input Parameter: 1809 . draw - the drawing context 1810 1811 Level: intermediate 1812 1813 .keywords: TS, monitor, line graph, destroy 1814 1815 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1816 @*/ 1817 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 1818 { 1819 PetscDraw draw; 1820 PetscErrorCode ierr; 1821 1822 PetscFunctionBegin; 1823 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 1824 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 1825 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1826 PetscFunctionReturn(0); 1827 } 1828 1829 #undef __FUNCT__ 1830 #define __FUNCT__ "TSGetTime" 1831 /*@ 1832 TSGetTime - Gets the current time. 1833 1834 Not Collective 1835 1836 Input Parameter: 1837 . ts - the TS context obtained from TSCreate() 1838 1839 Output Parameter: 1840 . t - the current time 1841 1842 Level: beginner 1843 1844 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1845 1846 .keywords: TS, get, time 1847 @*/ 1848 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1849 { 1850 PetscFunctionBegin; 1851 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1852 PetscValidDoublePointer(t,2); 1853 *t = ts->ptime; 1854 PetscFunctionReturn(0); 1855 } 1856 1857 #undef __FUNCT__ 1858 #define __FUNCT__ "TSSetTime" 1859 /*@ 1860 TSSetTime - Allows one to reset the time. 1861 1862 Logically Collective on TS 1863 1864 Input Parameters: 1865 + ts - the TS context obtained from TSCreate() 1866 - time - the time 1867 1868 Level: intermediate 1869 1870 .seealso: TSGetTime(), TSSetDuration() 1871 1872 .keywords: TS, set, time 1873 @*/ 1874 PetscErrorCode TSSetTime(TS ts, PetscReal t) 1875 { 1876 PetscFunctionBegin; 1877 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1878 PetscValidLogicalCollectiveReal(ts,t,2); 1879 ts->ptime = t; 1880 PetscFunctionReturn(0); 1881 } 1882 1883 #undef __FUNCT__ 1884 #define __FUNCT__ "TSSetOptionsPrefix" 1885 /*@C 1886 TSSetOptionsPrefix - Sets the prefix used for searching for all 1887 TS options in the database. 1888 1889 Logically Collective on TS 1890 1891 Input Parameter: 1892 + ts - The TS context 1893 - prefix - The prefix to prepend to all option names 1894 1895 Notes: 1896 A hyphen (-) must NOT be given at the beginning of the prefix name. 1897 The first character of all runtime options is AUTOMATICALLY the 1898 hyphen. 1899 1900 Level: advanced 1901 1902 .keywords: TS, set, options, prefix, database 1903 1904 .seealso: TSSetFromOptions() 1905 1906 @*/ 1907 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1908 { 1909 PetscErrorCode ierr; 1910 SNES snes; 1911 1912 PetscFunctionBegin; 1913 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1914 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1915 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1916 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1917 PetscFunctionReturn(0); 1918 } 1919 1920 1921 #undef __FUNCT__ 1922 #define __FUNCT__ "TSAppendOptionsPrefix" 1923 /*@C 1924 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1925 TS options in the database. 1926 1927 Logically Collective on TS 1928 1929 Input Parameter: 1930 + ts - The TS context 1931 - prefix - The prefix to prepend to all option names 1932 1933 Notes: 1934 A hyphen (-) must NOT be given at the beginning of the prefix name. 1935 The first character of all runtime options is AUTOMATICALLY the 1936 hyphen. 1937 1938 Level: advanced 1939 1940 .keywords: TS, append, options, prefix, database 1941 1942 .seealso: TSGetOptionsPrefix() 1943 1944 @*/ 1945 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1946 { 1947 PetscErrorCode ierr; 1948 SNES snes; 1949 1950 PetscFunctionBegin; 1951 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1952 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1953 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1954 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1955 PetscFunctionReturn(0); 1956 } 1957 1958 #undef __FUNCT__ 1959 #define __FUNCT__ "TSGetOptionsPrefix" 1960 /*@C 1961 TSGetOptionsPrefix - Sets the prefix used for searching for all 1962 TS options in the database. 1963 1964 Not Collective 1965 1966 Input Parameter: 1967 . ts - The TS context 1968 1969 Output Parameter: 1970 . prefix - A pointer to the prefix string used 1971 1972 Notes: On the fortran side, the user should pass in a string 'prifix' of 1973 sufficient length to hold the prefix. 1974 1975 Level: intermediate 1976 1977 .keywords: TS, get, options, prefix, database 1978 1979 .seealso: TSAppendOptionsPrefix() 1980 @*/ 1981 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 1982 { 1983 PetscErrorCode ierr; 1984 1985 PetscFunctionBegin; 1986 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1987 PetscValidPointer(prefix,2); 1988 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1989 PetscFunctionReturn(0); 1990 } 1991 1992 #undef __FUNCT__ 1993 #define __FUNCT__ "TSGetRHSJacobian" 1994 /*@C 1995 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 1996 1997 Not Collective, but parallel objects are returned if TS is parallel 1998 1999 Input Parameter: 2000 . ts - The TS context obtained from TSCreate() 2001 2002 Output Parameters: 2003 + J - The Jacobian J of F, where U_t = F(U,t) 2004 . M - The preconditioner matrix, usually the same as J 2005 . func - Function to compute the Jacobian of the RHS 2006 - ctx - User-defined context for Jacobian evaluation routine 2007 2008 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2009 2010 Level: intermediate 2011 2012 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2013 2014 .keywords: TS, timestep, get, matrix, Jacobian 2015 @*/ 2016 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2017 { 2018 PetscErrorCode ierr; 2019 SNES snes; 2020 2021 PetscFunctionBegin; 2022 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2023 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2024 if (func) *func = ts->ops->rhsjacobian; 2025 if (ctx) *ctx = ts->jacP; 2026 PetscFunctionReturn(0); 2027 } 2028 2029 #undef __FUNCT__ 2030 #define __FUNCT__ "TSGetIJacobian" 2031 /*@C 2032 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2033 2034 Not Collective, but parallel objects are returned if TS is parallel 2035 2036 Input Parameter: 2037 . ts - The TS context obtained from TSCreate() 2038 2039 Output Parameters: 2040 + A - The Jacobian of F(t,U,U_t) 2041 . B - The preconditioner matrix, often the same as A 2042 . f - The function to compute the matrices 2043 - ctx - User-defined context for Jacobian evaluation routine 2044 2045 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2046 2047 Level: advanced 2048 2049 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2050 2051 .keywords: TS, timestep, get, matrix, Jacobian 2052 @*/ 2053 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2054 { 2055 PetscErrorCode ierr; 2056 SNES snes; 2057 2058 PetscFunctionBegin; 2059 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2060 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2061 if (f) *f = ts->ops->ijacobian; 2062 if (ctx) *ctx = ts->jacP; 2063 PetscFunctionReturn(0); 2064 } 2065 2066 #undef __FUNCT__ 2067 #define __FUNCT__ "TSMonitorSolution" 2068 /*@C 2069 TSMonitorSolution - Monitors progress of the TS solvers by calling 2070 VecView() for the solution at each timestep 2071 2072 Collective on TS 2073 2074 Input Parameters: 2075 + ts - the TS context 2076 . step - current time-step 2077 . ptime - current time 2078 - dummy - either a viewer or PETSC_NULL 2079 2080 Level: intermediate 2081 2082 .keywords: TS, vector, monitor, view 2083 2084 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2085 @*/ 2086 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2087 { 2088 PetscErrorCode ierr; 2089 PetscViewer viewer = (PetscViewer) dummy; 2090 2091 PetscFunctionBegin; 2092 if (!dummy) { 2093 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2094 } 2095 ierr = VecView(x,viewer);CHKERRQ(ierr); 2096 PetscFunctionReturn(0); 2097 } 2098 2099 2100 #undef __FUNCT__ 2101 #define __FUNCT__ "TSSetDM" 2102 /*@ 2103 TSSetDM - Sets the DM that may be used by some preconditioners 2104 2105 Logically Collective on TS and DM 2106 2107 Input Parameters: 2108 + ts - the preconditioner context 2109 - dm - the dm 2110 2111 Level: intermediate 2112 2113 2114 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2115 @*/ 2116 PetscErrorCode TSSetDM(TS ts,DM dm) 2117 { 2118 PetscErrorCode ierr; 2119 SNES snes; 2120 2121 PetscFunctionBegin; 2122 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2123 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2124 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2125 ts->dm = dm; 2126 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2127 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2128 PetscFunctionReturn(0); 2129 } 2130 2131 #undef __FUNCT__ 2132 #define __FUNCT__ "TSGetDM" 2133 /*@ 2134 TSGetDM - Gets the DM that may be used by some preconditioners 2135 2136 Not Collective 2137 2138 Input Parameter: 2139 . ts - the preconditioner context 2140 2141 Output Parameter: 2142 . dm - the dm 2143 2144 Level: intermediate 2145 2146 2147 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2148 @*/ 2149 PetscErrorCode TSGetDM(TS ts,DM *dm) 2150 { 2151 PetscFunctionBegin; 2152 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2153 *dm = ts->dm; 2154 PetscFunctionReturn(0); 2155 } 2156 2157 #undef __FUNCT__ 2158 #define __FUNCT__ "SNESTSFormFunction" 2159 /*@ 2160 SNESTSFormFunction - Function to evaluate nonlinear residual 2161 2162 Logically Collective on SNES 2163 2164 Input Parameter: 2165 + snes - nonlinear solver 2166 . X - the current state at which to evaluate the residual 2167 - ctx - user context, must be a TS 2168 2169 Output Parameter: 2170 . F - the nonlinear residual 2171 2172 Notes: 2173 This function is not normally called by users and is automatically registered with the SNES used by TS. 2174 It is most frequently passed to MatFDColoringSetFunction(). 2175 2176 Level: advanced 2177 2178 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2179 @*/ 2180 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2181 { 2182 TS ts = (TS)ctx; 2183 PetscErrorCode ierr; 2184 2185 PetscFunctionBegin; 2186 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2187 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2188 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2189 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2190 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2191 PetscFunctionReturn(0); 2192 } 2193 2194 #undef __FUNCT__ 2195 #define __FUNCT__ "SNESTSFormJacobian" 2196 /*@ 2197 SNESTSFormJacobian - Function to evaluate the Jacobian 2198 2199 Collective on SNES 2200 2201 Input Parameter: 2202 + snes - nonlinear solver 2203 . X - the current state at which to evaluate the residual 2204 - ctx - user context, must be a TS 2205 2206 Output Parameter: 2207 + A - the Jacobian 2208 . B - the preconditioning matrix (may be the same as A) 2209 - flag - indicates any structure change in the matrix 2210 2211 Notes: 2212 This function is not normally called by users and is automatically registered with the SNES used by TS. 2213 2214 Level: developer 2215 2216 .seealso: SNESSetJacobian() 2217 @*/ 2218 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2219 { 2220 TS ts = (TS)ctx; 2221 PetscErrorCode ierr; 2222 2223 PetscFunctionBegin; 2224 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2225 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2226 PetscValidPointer(A,3); 2227 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2228 PetscValidPointer(B,4); 2229 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2230 PetscValidPointer(flag,5); 2231 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2232 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2233 PetscFunctionReturn(0); 2234 } 2235 2236 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2237 #include <mex.h> 2238 2239 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2240 2241 #undef __FUNCT__ 2242 #define __FUNCT__ "TSComputeFunction_Matlab" 2243 /* 2244 TSComputeFunction_Matlab - Calls the function that has been set with 2245 TSSetFunctionMatlab(). 2246 2247 Collective on TS 2248 2249 Input Parameters: 2250 + snes - the TS context 2251 - x - input vector 2252 2253 Output Parameter: 2254 . y - function vector, as set by TSSetFunction() 2255 2256 Notes: 2257 TSComputeFunction() is typically used within nonlinear solvers 2258 implementations, so most users would not generally call this routine 2259 themselves. 2260 2261 Level: developer 2262 2263 .keywords: TS, nonlinear, compute, function 2264 2265 .seealso: TSSetFunction(), TSGetFunction() 2266 */ 2267 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2268 { 2269 PetscErrorCode ierr; 2270 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2271 int nlhs = 1,nrhs = 7; 2272 mxArray *plhs[1],*prhs[7]; 2273 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2274 2275 PetscFunctionBegin; 2276 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2277 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2278 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2279 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2280 PetscCheckSameComm(snes,1,x,3); 2281 PetscCheckSameComm(snes,1,y,5); 2282 2283 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2284 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2285 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2286 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2287 prhs[0] = mxCreateDoubleScalar((double)ls); 2288 prhs[1] = mxCreateDoubleScalar(time); 2289 prhs[2] = mxCreateDoubleScalar((double)lx); 2290 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2291 prhs[4] = mxCreateDoubleScalar((double)ly); 2292 prhs[5] = mxCreateString(sctx->funcname); 2293 prhs[6] = sctx->ctx; 2294 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2295 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2296 mxDestroyArray(prhs[0]); 2297 mxDestroyArray(prhs[1]); 2298 mxDestroyArray(prhs[2]); 2299 mxDestroyArray(prhs[3]); 2300 mxDestroyArray(prhs[4]); 2301 mxDestroyArray(prhs[5]); 2302 mxDestroyArray(plhs[0]); 2303 PetscFunctionReturn(0); 2304 } 2305 2306 2307 #undef __FUNCT__ 2308 #define __FUNCT__ "TSSetFunctionMatlab" 2309 /* 2310 TSSetFunctionMatlab - Sets the function evaluation routine and function 2311 vector for use by the TS routines in solving ODEs 2312 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2313 2314 Logically Collective on TS 2315 2316 Input Parameters: 2317 + ts - the TS context 2318 - func - function evaluation routine 2319 2320 Calling sequence of func: 2321 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2322 2323 Level: beginner 2324 2325 .keywords: TS, nonlinear, set, function 2326 2327 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2328 */ 2329 PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2330 { 2331 PetscErrorCode ierr; 2332 TSMatlabContext *sctx; 2333 2334 PetscFunctionBegin; 2335 /* currently sctx is memory bleed */ 2336 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2337 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2338 /* 2339 This should work, but it doesn't 2340 sctx->ctx = ctx; 2341 mexMakeArrayPersistent(sctx->ctx); 2342 */ 2343 sctx->ctx = mxDuplicateArray(ctx); 2344 ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2345 PetscFunctionReturn(0); 2346 } 2347 2348 #undef __FUNCT__ 2349 #define __FUNCT__ "TSComputeJacobian_Matlab" 2350 /* 2351 TSComputeJacobian_Matlab - Calls the function that has been set with 2352 TSSetJacobianMatlab(). 2353 2354 Collective on TS 2355 2356 Input Parameters: 2357 + snes - the TS context 2358 . x - input vector 2359 . A, B - the matrices 2360 - ctx - user context 2361 2362 Output Parameter: 2363 . flag - structure of the matrix 2364 2365 Level: developer 2366 2367 .keywords: TS, nonlinear, compute, function 2368 2369 .seealso: TSSetFunction(), TSGetFunction() 2370 @*/ 2371 PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2372 { 2373 PetscErrorCode ierr; 2374 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2375 int nlhs = 2,nrhs = 9; 2376 mxArray *plhs[2],*prhs[9]; 2377 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2378 2379 PetscFunctionBegin; 2380 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2381 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2382 2383 /* call Matlab function in ctx with arguments x and y */ 2384 2385 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2386 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2387 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2388 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2389 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2390 prhs[0] = mxCreateDoubleScalar((double)ls); 2391 prhs[1] = mxCreateDoubleScalar((double)time); 2392 prhs[2] = mxCreateDoubleScalar((double)lx); 2393 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2394 prhs[4] = mxCreateDoubleScalar((double)shift); 2395 prhs[5] = mxCreateDoubleScalar((double)lA); 2396 prhs[6] = mxCreateDoubleScalar((double)lB); 2397 prhs[7] = mxCreateString(sctx->funcname); 2398 prhs[8] = sctx->ctx; 2399 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2400 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2401 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2402 mxDestroyArray(prhs[0]); 2403 mxDestroyArray(prhs[1]); 2404 mxDestroyArray(prhs[2]); 2405 mxDestroyArray(prhs[3]); 2406 mxDestroyArray(prhs[4]); 2407 mxDestroyArray(prhs[5]); 2408 mxDestroyArray(prhs[6]); 2409 mxDestroyArray(prhs[7]); 2410 mxDestroyArray(plhs[0]); 2411 mxDestroyArray(plhs[1]); 2412 PetscFunctionReturn(0); 2413 } 2414 2415 2416 #undef __FUNCT__ 2417 #define __FUNCT__ "TSSetJacobianMatlab" 2418 /* 2419 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2420 vector for use by the TS routines in solving ODEs from MATLAB. Here the function is a string containing the name of a MATLAB function 2421 2422 Logically Collective on TS 2423 2424 Input Parameters: 2425 + snes - the TS context 2426 . A,B - Jacobian matrices 2427 . func - function evaluation routine 2428 - ctx - user context 2429 2430 Calling sequence of func: 2431 $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2432 2433 2434 Level: developer 2435 2436 .keywords: TS, nonlinear, set, function 2437 2438 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2439 */ 2440 PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2441 { 2442 PetscErrorCode ierr; 2443 TSMatlabContext *sctx; 2444 2445 PetscFunctionBegin; 2446 /* currently sctx is memory bleed */ 2447 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2448 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2449 /* 2450 This should work, but it doesn't 2451 sctx->ctx = ctx; 2452 mexMakeArrayPersistent(sctx->ctx); 2453 */ 2454 sctx->ctx = mxDuplicateArray(ctx); 2455 ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2456 PetscFunctionReturn(0); 2457 } 2458 2459 #undef __FUNCT__ 2460 #define __FUNCT__ "TSMonitor_Matlab" 2461 /* 2462 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2463 2464 Collective on TS 2465 2466 .seealso: TSSetFunction(), TSGetFunction() 2467 @*/ 2468 PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2469 { 2470 PetscErrorCode ierr; 2471 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2472 int nlhs = 1,nrhs = 6; 2473 mxArray *plhs[1],*prhs[6]; 2474 long long int lx = 0,ls = 0; 2475 2476 PetscFunctionBegin; 2477 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2478 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2479 2480 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2481 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2482 prhs[0] = mxCreateDoubleScalar((double)ls); 2483 prhs[1] = mxCreateDoubleScalar((double)it); 2484 prhs[2] = mxCreateDoubleScalar((double)time); 2485 prhs[3] = mxCreateDoubleScalar((double)lx); 2486 prhs[4] = mxCreateString(sctx->funcname); 2487 prhs[5] = sctx->ctx; 2488 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2489 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2490 mxDestroyArray(prhs[0]); 2491 mxDestroyArray(prhs[1]); 2492 mxDestroyArray(prhs[2]); 2493 mxDestroyArray(prhs[3]); 2494 mxDestroyArray(prhs[4]); 2495 mxDestroyArray(plhs[0]); 2496 PetscFunctionReturn(0); 2497 } 2498 2499 2500 #undef __FUNCT__ 2501 #define __FUNCT__ "TSMonitorSetMatlab" 2502 /* 2503 TSMonitorSetMatlab - Sets the monitor function from Matlab 2504 2505 Level: developer 2506 2507 .keywords: TS, nonlinear, set, function 2508 2509 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2510 */ 2511 PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2512 { 2513 PetscErrorCode ierr; 2514 TSMatlabContext *sctx; 2515 2516 PetscFunctionBegin; 2517 /* currently sctx is memory bleed */ 2518 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2519 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2520 /* 2521 This should work, but it doesn't 2522 sctx->ctx = ctx; 2523 mexMakeArrayPersistent(sctx->ctx); 2524 */ 2525 sctx->ctx = mxDuplicateArray(ctx); 2526 ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2527 PetscFunctionReturn(0); 2528 } 2529 2530 #endif 2531