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