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