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