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