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