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