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