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 ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr); 93 ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections","",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr); 94 ierr = PetscOptionsBool("-ts_error_if_step_failed","Error if no step succeeds","",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr); 95 96 /* Monitor options */ 97 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 98 if (flg) { 99 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 100 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 101 } 102 opt = PETSC_FALSE; 103 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 104 if (opt) { 105 ierr = TSMonitorSet(ts,TSMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 106 } 107 opt = PETSC_FALSE; 108 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 109 if (opt) { 110 ierr = TSMonitorSet(ts,TSMonitorSolution,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 111 } 112 113 /* Handle TS type options */ 114 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 115 116 /* Handle specific TS options */ 117 if (ts->ops->setfromoptions) { 118 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 119 } 120 121 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 122 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 123 ierr = PetscOptionsEnd();CHKERRQ(ierr); 124 125 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 126 /* Handle subobject options */ 127 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 128 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 129 PetscFunctionReturn(0); 130 } 131 132 #undef __FUNCT__ 133 #define __FUNCT__ "TSViewFromOptions" 134 /*@ 135 TSViewFromOptions - This function visualizes the ts based upon user options. 136 137 Collective on TS 138 139 Input Parameter: 140 . ts - The ts 141 142 Level: intermediate 143 144 .keywords: TS, view, options, database 145 .seealso: TSSetFromOptions(), TSView() 146 @*/ 147 PetscErrorCode TSViewFromOptions(TS ts,const char title[]) 148 { 149 PetscViewer viewer; 150 PetscDraw draw; 151 PetscBool opt = PETSC_FALSE; 152 char fileName[PETSC_MAX_PATH_LEN]; 153 PetscErrorCode ierr; 154 155 PetscFunctionBegin; 156 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix, "-ts_view", fileName, PETSC_MAX_PATH_LEN, &opt);CHKERRQ(ierr); 157 if (opt && !PetscPreLoadingOn) { 158 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,fileName,&viewer);CHKERRQ(ierr); 159 ierr = TSView(ts, viewer);CHKERRQ(ierr); 160 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 161 } 162 opt = PETSC_FALSE; 163 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix, "-ts_view_draw", &opt,PETSC_NULL);CHKERRQ(ierr); 164 if (opt) { 165 ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm, 0, 0, 0, 0, 300, 300, &viewer);CHKERRQ(ierr); 166 ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 167 if (title) { 168 ierr = PetscDrawSetTitle(draw, (char *)title);CHKERRQ(ierr); 169 } else { 170 ierr = PetscObjectName((PetscObject)ts);CHKERRQ(ierr); 171 ierr = PetscDrawSetTitle(draw, ((PetscObject)ts)->name);CHKERRQ(ierr); 172 } 173 ierr = TSView(ts, viewer);CHKERRQ(ierr); 174 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 175 ierr = PetscDrawPause(draw);CHKERRQ(ierr); 176 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 177 } 178 PetscFunctionReturn(0); 179 } 180 181 #undef __FUNCT__ 182 #define __FUNCT__ "TSComputeRHSJacobian" 183 /*@ 184 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 185 set with TSSetRHSJacobian(). 186 187 Collective on TS and Vec 188 189 Input Parameters: 190 + ts - the TS context 191 . t - current timestep 192 - x - input vector 193 194 Output Parameters: 195 + A - Jacobian matrix 196 . B - optional preconditioning matrix 197 - flag - flag indicating matrix structure 198 199 Notes: 200 Most users should not need to explicitly call this routine, as it 201 is used internally within the nonlinear solvers. 202 203 See KSPSetOperators() for important information about setting the 204 flag parameter. 205 206 Level: developer 207 208 .keywords: SNES, compute, Jacobian, matrix 209 210 .seealso: TSSetRHSJacobian(), KSPSetOperators() 211 @*/ 212 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 213 { 214 PetscErrorCode ierr; 215 216 PetscFunctionBegin; 217 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 218 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 219 PetscCheckSameComm(ts,1,X,3); 220 if (ts->ops->rhsjacobian) { 221 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 222 *flg = DIFFERENT_NONZERO_PATTERN; 223 PetscStackPush("TS user Jacobian function"); 224 ierr = (*ts->ops->rhsjacobian)(ts,t,X,A,B,flg,ts->jacP);CHKERRQ(ierr); 225 PetscStackPop; 226 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 227 /* make sure user returned a correct Jacobian and preconditioner */ 228 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 229 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 230 } else { 231 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 232 if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 233 *flg = SAME_NONZERO_PATTERN; 234 } 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "TSComputeRHSFunction" 240 /*@ 241 TSComputeRHSFunction - Evaluates the right-hand-side function. 242 243 Collective on TS and Vec 244 245 Input Parameters: 246 + ts - the TS context 247 . t - current time 248 - x - state vector 249 250 Output Parameter: 251 . y - right hand side 252 253 Note: 254 Most users should not need to explicitly call this routine, as it 255 is used internally within the nonlinear solvers. 256 257 Level: developer 258 259 .keywords: TS, compute 260 261 .seealso: TSSetRHSFunction(), TSComputeIFunction() 262 @*/ 263 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 264 { 265 PetscErrorCode ierr; 266 267 PetscFunctionBegin; 268 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 269 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 270 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 271 272 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 273 if (ts->ops->rhsfunction) { 274 PetscStackPush("TS user right-hand-side function"); 275 ierr = (*ts->ops->rhsfunction)(ts,t,x,y,ts->funP);CHKERRQ(ierr); 276 PetscStackPop; 277 } else { 278 ierr = VecZeroEntries(y);CHKERRQ(ierr); 279 } 280 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNCT__ 285 #define __FUNCT__ "TSGetRHSVec_Private" 286 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 287 { 288 PetscErrorCode ierr; 289 290 PetscFunctionBegin; 291 if (!ts->Frhs) { 292 ierr = VecDuplicate(ts->vec_sol,&ts->Frhs);CHKERRQ(ierr); 293 } 294 *Frhs = ts->Frhs; 295 PetscFunctionReturn(0); 296 } 297 298 #undef __FUNCT__ 299 #define __FUNCT__ "TSGetRHSMats_Private" 300 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 301 { 302 PetscErrorCode ierr; 303 Mat A,B; 304 305 PetscFunctionBegin; 306 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 307 if (Arhs) { 308 if (!ts->Arhs) { 309 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 310 } 311 *Arhs = ts->Arhs; 312 } 313 if (Brhs) { 314 if (!ts->Brhs) { 315 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 316 } 317 *Brhs = ts->Brhs; 318 } 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNCT__ 323 #define __FUNCT__ "TSComputeIFunction" 324 /*@ 325 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 326 327 Collective on TS and Vec 328 329 Input Parameters: 330 + ts - the TS context 331 . t - current time 332 . X - state vector 333 . Xdot - time derivative of state vector 334 - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 335 336 Output Parameter: 337 . Y - right hand side 338 339 Note: 340 Most users should not need to explicitly call this routine, as it 341 is used internally within the nonlinear solvers. 342 343 If the user did did not write their equations in implicit form, this 344 function recasts them in implicit form. 345 346 Level: developer 347 348 .keywords: TS, compute 349 350 .seealso: TSSetIFunction(), TSComputeRHSFunction() 351 @*/ 352 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 353 { 354 PetscErrorCode ierr; 355 356 PetscFunctionBegin; 357 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 358 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 359 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 360 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 361 362 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 363 if (ts->ops->ifunction) { 364 PetscStackPush("TS user implicit function"); 365 ierr = (*ts->ops->ifunction)(ts,t,X,Xdot,Y,ts->funP);CHKERRQ(ierr); 366 PetscStackPop; 367 } 368 if (imex) { 369 if (!ts->ops->ifunction) {ierr = VecCopy(Xdot,Y);CHKERRQ(ierr);} 370 } else { 371 if (!ts->ops->ifunction) { 372 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 373 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 374 } else { 375 Vec Frhs; 376 ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 377 ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 378 ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 379 } 380 } 381 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "TSComputeIJacobian" 387 /*@ 388 TSComputeIJacobian - Evaluates the Jacobian of the DAE 389 390 Collective on TS and Vec 391 392 Input 393 Input Parameters: 394 + ts - the TS context 395 . t - current timestep 396 . X - state vector 397 . Xdot - time derivative of state vector 398 . shift - shift to apply, see note below 399 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 400 401 Output Parameters: 402 + A - Jacobian matrix 403 . B - optional preconditioning matrix 404 - flag - flag indicating matrix structure 405 406 Notes: 407 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 408 409 dF/dX + shift*dF/dXdot 410 411 Most users should not need to explicitly call this routine, as it 412 is used internally within the nonlinear solvers. 413 414 Level: developer 415 416 .keywords: TS, compute, Jacobian, matrix 417 418 .seealso: TSSetIJacobian() 419 @*/ 420 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 421 { 422 PetscErrorCode ierr; 423 424 PetscFunctionBegin; 425 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 426 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 427 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 428 PetscValidPointer(A,6); 429 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 430 PetscValidPointer(B,7); 431 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 432 PetscValidPointer(flg,8); 433 434 *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 435 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 436 if (ts->ops->ijacobian) { 437 PetscStackPush("TS user implicit Jacobian"); 438 ierr = (*ts->ops->ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ts->jacP);CHKERRQ(ierr); 439 PetscStackPop; 440 /* make sure user returned a correct Jacobian and preconditioner */ 441 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 442 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 443 } 444 if (imex) { 445 if (!ts->ops->ijacobian) { /* system was written as Xdot = F(t,X) */ 446 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 447 ierr = MatShift(*A,1.0);CHKERRQ(ierr); 448 if (*A != *B) { 449 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 450 ierr = MatShift(*B,shift);CHKERRQ(ierr); 451 } 452 *flg = SAME_PRECONDITIONER; 453 } 454 } else { 455 if (!ts->ops->ijacobian) { 456 ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 457 ierr = MatScale(*A,-1);CHKERRQ(ierr); 458 ierr = MatShift(*A,shift);CHKERRQ(ierr); 459 if (*A != *B) { 460 ierr = MatScale(*B,-1);CHKERRQ(ierr); 461 ierr = MatShift(*B,shift);CHKERRQ(ierr); 462 } 463 } else if (ts->ops->rhsjacobian) { 464 Mat Arhs,Brhs; 465 MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 466 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 467 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 468 axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 469 ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 470 if (*A != *B) { 471 ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 472 } 473 *flg = PetscMin(*flg,flg2); 474 } 475 } 476 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "TSSetRHSFunction" 482 /*@C 483 TSSetRHSFunction - Sets the routine for evaluating the function, 484 F(t,u), where U_t = F(t,u). 485 486 Logically Collective on TS 487 488 Input Parameters: 489 + ts - the TS context obtained from TSCreate() 490 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 491 . f - routine for evaluating the right-hand-side function 492 - ctx - [optional] user-defined context for private data for the 493 function evaluation routine (may be PETSC_NULL) 494 495 Calling sequence of func: 496 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 497 498 + t - current timestep 499 . u - input vector 500 . F - function vector 501 - ctx - [optional] user-defined function context 502 503 Important: 504 The user MUST call either this routine or TSSetMatrices(). 505 506 Level: beginner 507 508 .keywords: TS, timestep, set, right-hand-side, function 509 510 .seealso: TSSetMatrices() 511 @*/ 512 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 513 { 514 PetscErrorCode ierr; 515 SNES snes; 516 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 519 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 520 if (f) ts->ops->rhsfunction = f; 521 if (ctx) ts->funP = ctx; 522 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 523 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 524 PetscFunctionReturn(0); 525 } 526 527 #undef __FUNCT__ 528 #define __FUNCT__ "TSSetRHSJacobian" 529 /*@C 530 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 531 where U_t = F(U,t), as well as the location to store the matrix. 532 Use TSSetMatrices() for linear problems. 533 534 Logically Collective on TS 535 536 Input Parameters: 537 + ts - the TS context obtained from TSCreate() 538 . A - Jacobian matrix 539 . B - preconditioner matrix (usually same as A) 540 . f - the Jacobian evaluation routine 541 - ctx - [optional] user-defined context for private data for the 542 Jacobian evaluation routine (may be PETSC_NULL) 543 544 Calling sequence of func: 545 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 546 547 + t - current timestep 548 . u - input vector 549 . A - matrix A, where U_t = A(t)u 550 . B - preconditioner matrix, usually the same as A 551 . flag - flag indicating information about the preconditioner matrix 552 structure (same as flag in KSPSetOperators()) 553 - ctx - [optional] user-defined context for matrix evaluation routine 554 555 Notes: 556 See KSPSetOperators() for important information about setting the flag 557 output parameter in the routine func(). Be sure to read this information! 558 559 The routine func() takes Mat * as the matrix arguments rather than Mat. 560 This allows the matrix evaluation routine to replace A and/or B with a 561 completely new matrix structure (not just different matrix elements) 562 when appropriate, for instance, if the nonzero structure is changing 563 throughout the global iterations. 564 565 Level: beginner 566 567 .keywords: TS, timestep, set, right-hand-side, Jacobian 568 569 .seealso: TSDefaultComputeJacobianColor(), 570 SNESDefaultComputeJacobianColor(), TSSetRHSFunction(), TSSetMatrices() 571 572 @*/ 573 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 574 { 575 PetscErrorCode ierr; 576 SNES snes; 577 578 PetscFunctionBegin; 579 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 580 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 581 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 582 if (A) PetscCheckSameComm(ts,1,A,2); 583 if (B) PetscCheckSameComm(ts,1,B,3); 584 585 if (f) ts->ops->rhsjacobian = f; 586 if (ctx) ts->jacP = ctx; 587 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 588 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 589 PetscFunctionReturn(0); 590 } 591 592 593 #undef __FUNCT__ 594 #define __FUNCT__ "TSSetIFunction" 595 /*@C 596 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 597 598 Logically Collective on TS 599 600 Input Parameters: 601 + ts - the TS context obtained from TSCreate() 602 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 603 . f - the function evaluation routine 604 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 605 606 Calling sequence of f: 607 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 608 609 + t - time at step/stage being solved 610 . u - state vector 611 . u_t - time derivative of state vector 612 . F - function vector 613 - ctx - [optional] user-defined context for matrix evaluation routine 614 615 Important: 616 The user MUST call either this routine, TSSetRHSFunction(), or TSSetMatrices(). This routine must be used when not solving an ODE. 617 618 Level: beginner 619 620 .keywords: TS, timestep, set, DAE, Jacobian 621 622 .seealso: TSSetMatrices(), TSSetRHSFunction(), TSSetIJacobian() 623 @*/ 624 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 625 { 626 PetscErrorCode ierr; 627 SNES snes; 628 629 PetscFunctionBegin; 630 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 631 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 632 if (f) ts->ops->ifunction = f; 633 if (ctx) ts->funP = ctx; 634 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 635 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 636 PetscFunctionReturn(0); 637 } 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "TSGetIFunction" 641 /*@C 642 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 643 644 Not Collective 645 646 Input Parameter: 647 . ts - the TS context 648 649 Output Parameter: 650 + r - vector to hold residual (or PETSC_NULL) 651 . func - the function to compute residual (or PETSC_NULL) 652 - ctx - the function context (or PETSC_NULL) 653 654 Level: advanced 655 656 .keywords: TS, nonlinear, get, function 657 658 .seealso: TSSetIFunction(), SNESGetFunction() 659 @*/ 660 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 661 { 662 PetscErrorCode ierr; 663 SNES snes; 664 665 PetscFunctionBegin; 666 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 667 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 668 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 669 if (func) *func = ts->ops->ifunction; 670 if (ctx) *ctx = ts->funP; 671 PetscFunctionReturn(0); 672 } 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "TSGetRHSFunction" 676 /*@C 677 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 678 679 Not Collective 680 681 Input Parameter: 682 . ts - the TS context 683 684 Output Parameter: 685 + r - vector to hold computed right hand side (or PETSC_NULL) 686 . func - the function to compute right hand side (or PETSC_NULL) 687 - ctx - the function context (or PETSC_NULL) 688 689 Level: advanced 690 691 .keywords: TS, nonlinear, get, function 692 693 .seealso: TSSetRhsfunction(), SNESGetFunction() 694 @*/ 695 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 696 { 697 PetscErrorCode ierr; 698 SNES snes; 699 700 PetscFunctionBegin; 701 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 702 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 703 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 704 if (func) *func = ts->ops->rhsfunction; 705 if (ctx) *ctx = ts->funP; 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "TSSetIJacobian" 711 /*@C 712 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 713 you provided with TSSetIFunction(). 714 715 Logically Collective on TS 716 717 Input Parameters: 718 + ts - the TS context obtained from TSCreate() 719 . A - Jacobian matrix 720 . B - preconditioning matrix for A (may be same as A) 721 . f - the Jacobian evaluation routine 722 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 723 724 Calling sequence of f: 725 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 726 727 + t - time at step/stage being solved 728 . U - state vector 729 . U_t - time derivative of state vector 730 . a - shift 731 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 732 . B - preconditioning matrix for A, may be same as A 733 . flag - flag indicating information about the preconditioner matrix 734 structure (same as flag in KSPSetOperators()) 735 - ctx - [optional] user-defined context for matrix evaluation routine 736 737 Notes: 738 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 739 740 The matrix dF/dU + a*dF/dU_t you provide turns out to be 741 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 742 The time integrator internally approximates U_t by W+a*U where the positive "shift" 743 a and vector W depend on the integration method, step size, and past states. For example with 744 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 745 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 746 747 Level: beginner 748 749 .keywords: TS, timestep, DAE, Jacobian 750 751 .seealso: TSSetIFunction(), TSSetRHSJacobian() 752 753 @*/ 754 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 755 { 756 PetscErrorCode ierr; 757 SNES snes; 758 759 PetscFunctionBegin; 760 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 761 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 762 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 763 if (A) PetscCheckSameComm(ts,1,A,2); 764 if (B) PetscCheckSameComm(ts,1,B,3); 765 if (f) ts->ops->ijacobian = f; 766 if (ctx) ts->jacP = ctx; 767 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 768 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 769 PetscFunctionReturn(0); 770 } 771 772 #undef __FUNCT__ 773 #define __FUNCT__ "TSView" 774 /*@C 775 TSView - Prints the TS data structure. 776 777 Collective on TS 778 779 Input Parameters: 780 + ts - the TS context obtained from TSCreate() 781 - viewer - visualization context 782 783 Options Database Key: 784 . -ts_view - calls TSView() at end of TSStep() 785 786 Notes: 787 The available visualization contexts include 788 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 789 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 790 output where only the first processor opens 791 the file. All other processors send their 792 data to the first processor to print. 793 794 The user can open an alternative visualization context with 795 PetscViewerASCIIOpen() - output to a specified file. 796 797 Level: beginner 798 799 .keywords: TS, timestep, view 800 801 .seealso: PetscViewerASCIIOpen() 802 @*/ 803 PetscErrorCode TSView(TS ts,PetscViewer viewer) 804 { 805 PetscErrorCode ierr; 806 const TSType type; 807 PetscBool iascii,isstring,isundials; 808 809 PetscFunctionBegin; 810 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 811 if (!viewer) { 812 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 813 } 814 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 815 PetscCheckSameComm(ts,1,viewer,2); 816 817 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 818 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 819 if (iascii) { 820 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 821 if (ts->ops->view) { 822 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 823 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 824 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 825 } 826 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 827 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 828 if (ts->problem_type == TS_NONLINEAR) { 829 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->nonlinear_its);CHKERRQ(ierr); 830 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solve failures=%D\n",ts->max_snes_failures);CHKERRQ(ierr); 831 } 832 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->linear_its);CHKERRQ(ierr); 833 ierr = PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr); 834 } else if (isstring) { 835 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 836 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 837 } 838 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 839 ierr = PetscTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 840 if (!isundials && ts->snes) {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);} 841 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 842 PetscFunctionReturn(0); 843 } 844 845 846 #undef __FUNCT__ 847 #define __FUNCT__ "TSSetApplicationContext" 848 /*@ 849 TSSetApplicationContext - Sets an optional user-defined context for 850 the timesteppers. 851 852 Logically Collective on TS 853 854 Input Parameters: 855 + ts - the TS context obtained from TSCreate() 856 - usrP - optional user context 857 858 Level: intermediate 859 860 .keywords: TS, timestep, set, application, context 861 862 .seealso: TSGetApplicationContext() 863 @*/ 864 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 865 { 866 PetscFunctionBegin; 867 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 868 ts->user = usrP; 869 PetscFunctionReturn(0); 870 } 871 872 #undef __FUNCT__ 873 #define __FUNCT__ "TSGetApplicationContext" 874 /*@ 875 TSGetApplicationContext - Gets the user-defined context for the 876 timestepper. 877 878 Not Collective 879 880 Input Parameter: 881 . ts - the TS context obtained from TSCreate() 882 883 Output Parameter: 884 . usrP - user context 885 886 Level: intermediate 887 888 .keywords: TS, timestep, get, application, context 889 890 .seealso: TSSetApplicationContext() 891 @*/ 892 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 893 { 894 PetscFunctionBegin; 895 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 896 *(void**)usrP = ts->user; 897 PetscFunctionReturn(0); 898 } 899 900 #undef __FUNCT__ 901 #define __FUNCT__ "TSGetTimeStepNumber" 902 /*@ 903 TSGetTimeStepNumber - Gets the current number of timesteps. 904 905 Not Collective 906 907 Input Parameter: 908 . ts - the TS context obtained from TSCreate() 909 910 Output Parameter: 911 . iter - number steps so far 912 913 Level: intermediate 914 915 .keywords: TS, timestep, get, iteration, number 916 @*/ 917 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 918 { 919 PetscFunctionBegin; 920 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 921 PetscValidIntPointer(iter,2); 922 *iter = ts->steps; 923 PetscFunctionReturn(0); 924 } 925 926 #undef __FUNCT__ 927 #define __FUNCT__ "TSSetInitialTimeStep" 928 /*@ 929 TSSetInitialTimeStep - Sets the initial timestep to be used, 930 as well as the initial time. 931 932 Logically Collective on TS 933 934 Input Parameters: 935 + ts - the TS context obtained from TSCreate() 936 . initial_time - the initial time 937 - time_step - the size of the timestep 938 939 Level: intermediate 940 941 .seealso: TSSetTimeStep(), TSGetTimeStep() 942 943 .keywords: TS, set, initial, timestep 944 @*/ 945 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 946 { 947 PetscFunctionBegin; 948 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 949 ts->time_step = time_step; 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 /* reset time step and iteration counters */ 1704 ts->steps = 0; 1705 ts->linear_its = 0; 1706 ts->nonlinear_its = 0; 1707 ts->reason = TS_CONVERGED_ITERATING; 1708 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1709 1710 if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */ 1711 ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr); 1712 } else { 1713 /* steps the requested number of timesteps. */ 1714 for (i=0; !ts->reason; i++) { 1715 ierr = TSPreStep(ts);CHKERRQ(ierr); 1716 ierr = TSStep(ts);CHKERRQ(ierr); 1717 if (ts->reason < 0) { 1718 if (ts->errorifstepfailed) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed"); 1719 } else if (i >= ts->max_steps) { 1720 ts->reason = TS_CONVERGED_ITS; 1721 } else if (ts->ptime >= ts->max_time) { 1722 ts->reason = TS_CONVERGED_TIME; 1723 } 1724 ierr = TSPostStep(ts);CHKERRQ(ierr); 1725 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1726 } 1727 } 1728 if (!PetscPreLoadingOn) { 1729 ierr = TSViewFromOptions(ts,((PetscObject)ts)->name);CHKERRQ(ierr); 1730 } 1731 PetscFunctionReturn(0); 1732 } 1733 1734 #undef __FUNCT__ 1735 #define __FUNCT__ "TSMonitor" 1736 /* 1737 Runs the user provided monitor routines, if they exists. 1738 */ 1739 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 1740 { 1741 PetscErrorCode ierr; 1742 PetscInt i,n = ts->numbermonitors; 1743 1744 PetscFunctionBegin; 1745 for (i=0; i<n; i++) { 1746 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 1747 } 1748 PetscFunctionReturn(0); 1749 } 1750 1751 /* ------------------------------------------------------------------------*/ 1752 1753 #undef __FUNCT__ 1754 #define __FUNCT__ "TSMonitorLGCreate" 1755 /*@C 1756 TSMonitorLGCreate - Creates a line graph context for use with 1757 TS to monitor convergence of preconditioned residual norms. 1758 1759 Collective on TS 1760 1761 Input Parameters: 1762 + host - the X display to open, or null for the local machine 1763 . label - the title to put in the title bar 1764 . x, y - the screen coordinates of the upper left coordinate of the window 1765 - m, n - the screen width and height in pixels 1766 1767 Output Parameter: 1768 . draw - the drawing context 1769 1770 Options Database Key: 1771 . -ts_monitor_draw - automatically sets line graph monitor 1772 1773 Notes: 1774 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 1775 1776 Level: intermediate 1777 1778 .keywords: TS, monitor, line graph, residual, seealso 1779 1780 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 1781 1782 @*/ 1783 PetscErrorCode TSMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1784 { 1785 PetscDraw win; 1786 PetscErrorCode ierr; 1787 1788 PetscFunctionBegin; 1789 ierr = PetscDrawCreate(PETSC_COMM_SELF,host,label,x,y,m,n,&win);CHKERRQ(ierr); 1790 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 1791 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 1792 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 1793 1794 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 1795 PetscFunctionReturn(0); 1796 } 1797 1798 #undef __FUNCT__ 1799 #define __FUNCT__ "TSMonitorLG" 1800 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 1801 { 1802 PetscDrawLG lg = (PetscDrawLG) monctx; 1803 PetscReal x,y = ptime; 1804 PetscErrorCode ierr; 1805 1806 PetscFunctionBegin; 1807 if (!monctx) { 1808 MPI_Comm comm; 1809 PetscViewer viewer; 1810 1811 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 1812 viewer = PETSC_VIEWER_DRAW_(comm); 1813 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 1814 } 1815 1816 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 1817 x = (PetscReal)n; 1818 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 1819 if (n < 20 || (n % 5)) { 1820 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 1821 } 1822 PetscFunctionReturn(0); 1823 } 1824 1825 #undef __FUNCT__ 1826 #define __FUNCT__ "TSMonitorLGDestroy" 1827 /*@C 1828 TSMonitorLGDestroy - Destroys a line graph context that was created 1829 with TSMonitorLGCreate(). 1830 1831 Collective on PetscDrawLG 1832 1833 Input Parameter: 1834 . draw - the drawing context 1835 1836 Level: intermediate 1837 1838 .keywords: TS, monitor, line graph, destroy 1839 1840 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 1841 @*/ 1842 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 1843 { 1844 PetscDraw draw; 1845 PetscErrorCode ierr; 1846 1847 PetscFunctionBegin; 1848 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 1849 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 1850 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 1851 PetscFunctionReturn(0); 1852 } 1853 1854 #undef __FUNCT__ 1855 #define __FUNCT__ "TSGetTime" 1856 /*@ 1857 TSGetTime - Gets the current time. 1858 1859 Not Collective 1860 1861 Input Parameter: 1862 . ts - the TS context obtained from TSCreate() 1863 1864 Output Parameter: 1865 . t - the current time 1866 1867 Level: beginner 1868 1869 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1870 1871 .keywords: TS, get, time 1872 @*/ 1873 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 1874 { 1875 PetscFunctionBegin; 1876 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1877 PetscValidDoublePointer(t,2); 1878 *t = ts->ptime; 1879 PetscFunctionReturn(0); 1880 } 1881 1882 #undef __FUNCT__ 1883 #define __FUNCT__ "TSSetTime" 1884 /*@ 1885 TSSetTime - Allows one to reset the time. 1886 1887 Logically Collective on TS 1888 1889 Input Parameters: 1890 + ts - the TS context obtained from TSCreate() 1891 - time - the time 1892 1893 Level: intermediate 1894 1895 .seealso: TSGetTime(), TSSetDuration() 1896 1897 .keywords: TS, set, time 1898 @*/ 1899 PetscErrorCode TSSetTime(TS ts, PetscReal t) 1900 { 1901 PetscFunctionBegin; 1902 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1903 PetscValidLogicalCollectiveReal(ts,t,2); 1904 ts->ptime = t; 1905 PetscFunctionReturn(0); 1906 } 1907 1908 #undef __FUNCT__ 1909 #define __FUNCT__ "TSSetOptionsPrefix" 1910 /*@C 1911 TSSetOptionsPrefix - Sets the prefix used for searching for all 1912 TS options in the database. 1913 1914 Logically Collective on TS 1915 1916 Input Parameter: 1917 + ts - The TS context 1918 - prefix - The prefix to prepend to all option names 1919 1920 Notes: 1921 A hyphen (-) must NOT be given at the beginning of the prefix name. 1922 The first character of all runtime options is AUTOMATICALLY the 1923 hyphen. 1924 1925 Level: advanced 1926 1927 .keywords: TS, set, options, prefix, database 1928 1929 .seealso: TSSetFromOptions() 1930 1931 @*/ 1932 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 1933 { 1934 PetscErrorCode ierr; 1935 SNES snes; 1936 1937 PetscFunctionBegin; 1938 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1939 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1940 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1941 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1942 PetscFunctionReturn(0); 1943 } 1944 1945 1946 #undef __FUNCT__ 1947 #define __FUNCT__ "TSAppendOptionsPrefix" 1948 /*@C 1949 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 1950 TS options in the database. 1951 1952 Logically Collective on TS 1953 1954 Input Parameter: 1955 + ts - The TS context 1956 - prefix - The prefix to prepend to all option names 1957 1958 Notes: 1959 A hyphen (-) must NOT be given at the beginning of the prefix name. 1960 The first character of all runtime options is AUTOMATICALLY the 1961 hyphen. 1962 1963 Level: advanced 1964 1965 .keywords: TS, append, options, prefix, database 1966 1967 .seealso: TSGetOptionsPrefix() 1968 1969 @*/ 1970 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 1971 { 1972 PetscErrorCode ierr; 1973 SNES snes; 1974 1975 PetscFunctionBegin; 1976 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1977 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 1978 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1979 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 1980 PetscFunctionReturn(0); 1981 } 1982 1983 #undef __FUNCT__ 1984 #define __FUNCT__ "TSGetOptionsPrefix" 1985 /*@C 1986 TSGetOptionsPrefix - Sets the prefix used for searching for all 1987 TS options in the database. 1988 1989 Not Collective 1990 1991 Input Parameter: 1992 . ts - The TS context 1993 1994 Output Parameter: 1995 . prefix - A pointer to the prefix string used 1996 1997 Notes: On the fortran side, the user should pass in a string 'prifix' of 1998 sufficient length to hold the prefix. 1999 2000 Level: intermediate 2001 2002 .keywords: TS, get, options, prefix, database 2003 2004 .seealso: TSAppendOptionsPrefix() 2005 @*/ 2006 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2007 { 2008 PetscErrorCode ierr; 2009 2010 PetscFunctionBegin; 2011 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2012 PetscValidPointer(prefix,2); 2013 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2014 PetscFunctionReturn(0); 2015 } 2016 2017 #undef __FUNCT__ 2018 #define __FUNCT__ "TSGetRHSJacobian" 2019 /*@C 2020 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2021 2022 Not Collective, but parallel objects are returned if TS is parallel 2023 2024 Input Parameter: 2025 . ts - The TS context obtained from TSCreate() 2026 2027 Output Parameters: 2028 + J - The Jacobian J of F, where U_t = F(U,t) 2029 . M - The preconditioner matrix, usually the same as J 2030 . func - Function to compute the Jacobian of the RHS 2031 - ctx - User-defined context for Jacobian evaluation routine 2032 2033 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2034 2035 Level: intermediate 2036 2037 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2038 2039 .keywords: TS, timestep, get, matrix, Jacobian 2040 @*/ 2041 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2042 { 2043 PetscErrorCode ierr; 2044 SNES snes; 2045 2046 PetscFunctionBegin; 2047 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2048 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2049 if (func) *func = ts->ops->rhsjacobian; 2050 if (ctx) *ctx = ts->jacP; 2051 PetscFunctionReturn(0); 2052 } 2053 2054 #undef __FUNCT__ 2055 #define __FUNCT__ "TSGetIJacobian" 2056 /*@C 2057 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2058 2059 Not Collective, but parallel objects are returned if TS is parallel 2060 2061 Input Parameter: 2062 . ts - The TS context obtained from TSCreate() 2063 2064 Output Parameters: 2065 + A - The Jacobian of F(t,U,U_t) 2066 . B - The preconditioner matrix, often the same as A 2067 . f - The function to compute the matrices 2068 - ctx - User-defined context for Jacobian evaluation routine 2069 2070 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2071 2072 Level: advanced 2073 2074 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2075 2076 .keywords: TS, timestep, get, matrix, Jacobian 2077 @*/ 2078 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2079 { 2080 PetscErrorCode ierr; 2081 SNES snes; 2082 2083 PetscFunctionBegin; 2084 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2085 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2086 if (f) *f = ts->ops->ijacobian; 2087 if (ctx) *ctx = ts->jacP; 2088 PetscFunctionReturn(0); 2089 } 2090 2091 #undef __FUNCT__ 2092 #define __FUNCT__ "TSMonitorSolution" 2093 /*@C 2094 TSMonitorSolution - Monitors progress of the TS solvers by calling 2095 VecView() for the solution at each timestep 2096 2097 Collective on TS 2098 2099 Input Parameters: 2100 + ts - the TS context 2101 . step - current time-step 2102 . ptime - current time 2103 - dummy - either a viewer or PETSC_NULL 2104 2105 Level: intermediate 2106 2107 .keywords: TS, vector, monitor, view 2108 2109 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2110 @*/ 2111 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2112 { 2113 PetscErrorCode ierr; 2114 PetscViewer viewer = (PetscViewer) dummy; 2115 2116 PetscFunctionBegin; 2117 if (!dummy) { 2118 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2119 } 2120 ierr = VecView(x,viewer);CHKERRQ(ierr); 2121 PetscFunctionReturn(0); 2122 } 2123 2124 2125 #undef __FUNCT__ 2126 #define __FUNCT__ "TSSetDM" 2127 /*@ 2128 TSSetDM - Sets the DM that may be used by some preconditioners 2129 2130 Logically Collective on TS and DM 2131 2132 Input Parameters: 2133 + ts - the preconditioner context 2134 - dm - the dm 2135 2136 Level: intermediate 2137 2138 2139 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2140 @*/ 2141 PetscErrorCode TSSetDM(TS ts,DM dm) 2142 { 2143 PetscErrorCode ierr; 2144 SNES snes; 2145 2146 PetscFunctionBegin; 2147 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2148 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2149 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2150 ts->dm = dm; 2151 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2152 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2153 PetscFunctionReturn(0); 2154 } 2155 2156 #undef __FUNCT__ 2157 #define __FUNCT__ "TSGetDM" 2158 /*@ 2159 TSGetDM - Gets the DM that may be used by some preconditioners 2160 2161 Not Collective 2162 2163 Input Parameter: 2164 . ts - the preconditioner context 2165 2166 Output Parameter: 2167 . dm - the dm 2168 2169 Level: intermediate 2170 2171 2172 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2173 @*/ 2174 PetscErrorCode TSGetDM(TS ts,DM *dm) 2175 { 2176 PetscFunctionBegin; 2177 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2178 *dm = ts->dm; 2179 PetscFunctionReturn(0); 2180 } 2181 2182 #undef __FUNCT__ 2183 #define __FUNCT__ "SNESTSFormFunction" 2184 /*@ 2185 SNESTSFormFunction - Function to evaluate nonlinear residual 2186 2187 Logically Collective on SNES 2188 2189 Input Parameter: 2190 + snes - nonlinear solver 2191 . X - the current state at which to evaluate the residual 2192 - ctx - user context, must be a TS 2193 2194 Output Parameter: 2195 . F - the nonlinear residual 2196 2197 Notes: 2198 This function is not normally called by users and is automatically registered with the SNES used by TS. 2199 It is most frequently passed to MatFDColoringSetFunction(). 2200 2201 Level: advanced 2202 2203 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2204 @*/ 2205 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2206 { 2207 TS ts = (TS)ctx; 2208 PetscErrorCode ierr; 2209 2210 PetscFunctionBegin; 2211 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2212 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2213 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2214 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2215 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2216 PetscFunctionReturn(0); 2217 } 2218 2219 #undef __FUNCT__ 2220 #define __FUNCT__ "SNESTSFormJacobian" 2221 /*@ 2222 SNESTSFormJacobian - Function to evaluate the Jacobian 2223 2224 Collective on SNES 2225 2226 Input Parameter: 2227 + snes - nonlinear solver 2228 . X - the current state at which to evaluate the residual 2229 - ctx - user context, must be a TS 2230 2231 Output Parameter: 2232 + A - the Jacobian 2233 . B - the preconditioning matrix (may be the same as A) 2234 - flag - indicates any structure change in the matrix 2235 2236 Notes: 2237 This function is not normally called by users and is automatically registered with the SNES used by TS. 2238 2239 Level: developer 2240 2241 .seealso: SNESSetJacobian() 2242 @*/ 2243 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2244 { 2245 TS ts = (TS)ctx; 2246 PetscErrorCode ierr; 2247 2248 PetscFunctionBegin; 2249 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2250 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2251 PetscValidPointer(A,3); 2252 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2253 PetscValidPointer(B,4); 2254 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2255 PetscValidPointer(flag,5); 2256 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2257 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2258 PetscFunctionReturn(0); 2259 } 2260 2261 #if defined(PETSC_HAVE_MATLAB_ENGINE) 2262 #include <mex.h> 2263 2264 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 2265 2266 #undef __FUNCT__ 2267 #define __FUNCT__ "TSComputeFunction_Matlab" 2268 /* 2269 TSComputeFunction_Matlab - Calls the function that has been set with 2270 TSSetFunctionMatlab(). 2271 2272 Collective on TS 2273 2274 Input Parameters: 2275 + snes - the TS context 2276 - x - input vector 2277 2278 Output Parameter: 2279 . y - function vector, as set by TSSetFunction() 2280 2281 Notes: 2282 TSComputeFunction() is typically used within nonlinear solvers 2283 implementations, so most users would not generally call this routine 2284 themselves. 2285 2286 Level: developer 2287 2288 .keywords: TS, nonlinear, compute, function 2289 2290 .seealso: TSSetFunction(), TSGetFunction() 2291 */ 2292 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 2293 { 2294 PetscErrorCode ierr; 2295 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2296 int nlhs = 1,nrhs = 7; 2297 mxArray *plhs[1],*prhs[7]; 2298 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 2299 2300 PetscFunctionBegin; 2301 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2302 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2303 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 2304 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 2305 PetscCheckSameComm(snes,1,x,3); 2306 PetscCheckSameComm(snes,1,y,5); 2307 2308 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2309 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2310 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 2311 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 2312 prhs[0] = mxCreateDoubleScalar((double)ls); 2313 prhs[1] = mxCreateDoubleScalar(time); 2314 prhs[2] = mxCreateDoubleScalar((double)lx); 2315 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2316 prhs[4] = mxCreateDoubleScalar((double)ly); 2317 prhs[5] = mxCreateString(sctx->funcname); 2318 prhs[6] = sctx->ctx; 2319 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 2320 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2321 mxDestroyArray(prhs[0]); 2322 mxDestroyArray(prhs[1]); 2323 mxDestroyArray(prhs[2]); 2324 mxDestroyArray(prhs[3]); 2325 mxDestroyArray(prhs[4]); 2326 mxDestroyArray(prhs[5]); 2327 mxDestroyArray(plhs[0]); 2328 PetscFunctionReturn(0); 2329 } 2330 2331 2332 #undef __FUNCT__ 2333 #define __FUNCT__ "TSSetFunctionMatlab" 2334 /* 2335 TSSetFunctionMatlab - Sets the function evaluation routine and function 2336 vector for use by the TS routines in solving ODEs 2337 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 2338 2339 Logically Collective on TS 2340 2341 Input Parameters: 2342 + ts - the TS context 2343 - func - function evaluation routine 2344 2345 Calling sequence of func: 2346 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 2347 2348 Level: beginner 2349 2350 .keywords: TS, nonlinear, set, function 2351 2352 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2353 */ 2354 PetscErrorCode TSSetFunctionMatlab(TS snes,const char *func,mxArray *ctx) 2355 { 2356 PetscErrorCode ierr; 2357 TSMatlabContext *sctx; 2358 2359 PetscFunctionBegin; 2360 /* currently sctx is memory bleed */ 2361 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2362 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2363 /* 2364 This should work, but it doesn't 2365 sctx->ctx = ctx; 2366 mexMakeArrayPersistent(sctx->ctx); 2367 */ 2368 sctx->ctx = mxDuplicateArray(ctx); 2369 ierr = TSSetIFunction(snes,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 2370 PetscFunctionReturn(0); 2371 } 2372 2373 #undef __FUNCT__ 2374 #define __FUNCT__ "TSComputeJacobian_Matlab" 2375 /* 2376 TSComputeJacobian_Matlab - Calls the function that has been set with 2377 TSSetJacobianMatlab(). 2378 2379 Collective on TS 2380 2381 Input Parameters: 2382 + snes - the TS context 2383 . x - input vector 2384 . A, B - the matrices 2385 - ctx - user context 2386 2387 Output Parameter: 2388 . flag - structure of the matrix 2389 2390 Level: developer 2391 2392 .keywords: TS, nonlinear, compute, function 2393 2394 .seealso: TSSetFunction(), TSGetFunction() 2395 @*/ 2396 PetscErrorCode TSComputeJacobian_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 2397 { 2398 PetscErrorCode ierr; 2399 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2400 int nlhs = 2,nrhs = 9; 2401 mxArray *plhs[2],*prhs[9]; 2402 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 2403 2404 PetscFunctionBegin; 2405 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2406 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 2407 2408 /* call Matlab function in ctx with arguments x and y */ 2409 2410 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2411 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2412 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 2413 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 2414 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 2415 prhs[0] = mxCreateDoubleScalar((double)ls); 2416 prhs[1] = mxCreateDoubleScalar((double)time); 2417 prhs[2] = mxCreateDoubleScalar((double)lx); 2418 prhs[3] = mxCreateDoubleScalar((double)lxdot); 2419 prhs[4] = mxCreateDoubleScalar((double)shift); 2420 prhs[5] = mxCreateDoubleScalar((double)lA); 2421 prhs[6] = mxCreateDoubleScalar((double)lB); 2422 prhs[7] = mxCreateString(sctx->funcname); 2423 prhs[8] = sctx->ctx; 2424 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 2425 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2426 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 2427 mxDestroyArray(prhs[0]); 2428 mxDestroyArray(prhs[1]); 2429 mxDestroyArray(prhs[2]); 2430 mxDestroyArray(prhs[3]); 2431 mxDestroyArray(prhs[4]); 2432 mxDestroyArray(prhs[5]); 2433 mxDestroyArray(prhs[6]); 2434 mxDestroyArray(prhs[7]); 2435 mxDestroyArray(plhs[0]); 2436 mxDestroyArray(plhs[1]); 2437 PetscFunctionReturn(0); 2438 } 2439 2440 2441 #undef __FUNCT__ 2442 #define __FUNCT__ "TSSetJacobianMatlab" 2443 /* 2444 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 2445 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 2446 2447 Logically Collective on TS 2448 2449 Input Parameters: 2450 + snes - the TS context 2451 . A,B - Jacobian matrices 2452 . func - function evaluation routine 2453 - ctx - user context 2454 2455 Calling sequence of func: 2456 $ flag = func (TS snes,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 2457 2458 2459 Level: developer 2460 2461 .keywords: TS, nonlinear, set, function 2462 2463 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2464 */ 2465 PetscErrorCode TSSetJacobianMatlab(TS snes,Mat A,Mat B,const char *func,mxArray *ctx) 2466 { 2467 PetscErrorCode ierr; 2468 TSMatlabContext *sctx; 2469 2470 PetscFunctionBegin; 2471 /* currently sctx is memory bleed */ 2472 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2473 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2474 /* 2475 This should work, but it doesn't 2476 sctx->ctx = ctx; 2477 mexMakeArrayPersistent(sctx->ctx); 2478 */ 2479 sctx->ctx = mxDuplicateArray(ctx); 2480 ierr = TSSetIJacobian(snes,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 2481 PetscFunctionReturn(0); 2482 } 2483 2484 #undef __FUNCT__ 2485 #define __FUNCT__ "TSMonitor_Matlab" 2486 /* 2487 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 2488 2489 Collective on TS 2490 2491 .seealso: TSSetFunction(), TSGetFunction() 2492 @*/ 2493 PetscErrorCode TSMonitor_Matlab(TS snes,PetscInt it, PetscReal time,Vec x, void *ctx) 2494 { 2495 PetscErrorCode ierr; 2496 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 2497 int nlhs = 1,nrhs = 6; 2498 mxArray *plhs[1],*prhs[6]; 2499 long long int lx = 0,ls = 0; 2500 2501 PetscFunctionBegin; 2502 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 2503 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 2504 2505 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 2506 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 2507 prhs[0] = mxCreateDoubleScalar((double)ls); 2508 prhs[1] = mxCreateDoubleScalar((double)it); 2509 prhs[2] = mxCreateDoubleScalar((double)time); 2510 prhs[3] = mxCreateDoubleScalar((double)lx); 2511 prhs[4] = mxCreateString(sctx->funcname); 2512 prhs[5] = sctx->ctx; 2513 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 2514 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 2515 mxDestroyArray(prhs[0]); 2516 mxDestroyArray(prhs[1]); 2517 mxDestroyArray(prhs[2]); 2518 mxDestroyArray(prhs[3]); 2519 mxDestroyArray(prhs[4]); 2520 mxDestroyArray(plhs[0]); 2521 PetscFunctionReturn(0); 2522 } 2523 2524 2525 #undef __FUNCT__ 2526 #define __FUNCT__ "TSMonitorSetMatlab" 2527 /* 2528 TSMonitorSetMatlab - Sets the monitor function from Matlab 2529 2530 Level: developer 2531 2532 .keywords: TS, nonlinear, set, function 2533 2534 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 2535 */ 2536 PetscErrorCode TSMonitorSetMatlab(TS snes,const char *func,mxArray *ctx) 2537 { 2538 PetscErrorCode ierr; 2539 TSMatlabContext *sctx; 2540 2541 PetscFunctionBegin; 2542 /* currently sctx is memory bleed */ 2543 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 2544 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 2545 /* 2546 This should work, but it doesn't 2547 sctx->ctx = ctx; 2548 mexMakeArrayPersistent(sctx->ctx); 2549 */ 2550 sctx->ctx = mxDuplicateArray(ctx); 2551 ierr = TSMonitorSet(snes,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 2552 PetscFunctionReturn(0); 2553 } 2554 2555 #endif 2556