1 2 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 3 #include <petscdmshell.h> 4 5 /* Logging support */ 6 PetscClassId TS_CLASSID; 7 PetscLogEvent TS_Step, TS_PseudoComputeTimeStep, TS_FunctionEval, TS_JacobianEval; 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "TSSetTypeFromOptions" 11 /* 12 TSSetTypeFromOptions - Sets the type of ts from user options. 13 14 Collective on TS 15 16 Input Parameter: 17 . ts - The ts 18 19 Level: intermediate 20 21 .keywords: TS, set, options, database, type 22 .seealso: TSSetFromOptions(), TSSetType() 23 */ 24 static PetscErrorCode TSSetTypeFromOptions(TS ts) 25 { 26 PetscBool opt; 27 const char *defaultType; 28 char typeName[256]; 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 if (((PetscObject)ts)->type_name) { 33 defaultType = ((PetscObject)ts)->type_name; 34 } else { 35 defaultType = TSEULER; 36 } 37 38 if (!TSRegisterAllCalled) {ierr = TSRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 39 ierr = PetscOptionsList("-ts_type", "TS method"," TSSetType", TSList, defaultType, typeName, 256, &opt);CHKERRQ(ierr); 40 if (opt) { 41 ierr = TSSetType(ts, typeName);CHKERRQ(ierr); 42 } else { 43 ierr = TSSetType(ts, defaultType);CHKERRQ(ierr); 44 } 45 PetscFunctionReturn(0); 46 } 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "TSSetFromOptions" 50 /*@ 51 TSSetFromOptions - Sets various TS parameters from user options. 52 53 Collective on TS 54 55 Input Parameter: 56 . ts - the TS context obtained from TSCreate() 57 58 Options Database Keys: 59 + -ts_type <type> - TSEULER, TSBEULER, TSSUNDIALS, TSPSEUDO, TSCN, TSRK, TSTHETA, TSGL, TSSSP 60 . -ts_max_steps maxsteps - maximum number of time-steps to take 61 . -ts_final_time time - maximum time to compute to 62 . -ts_dt dt - initial time step 63 . -ts_monitor - print information at each timestep 64 - -ts_monitor_draw - plot information at each timestep 65 66 Level: beginner 67 68 .keywords: TS, timestep, set, options, database 69 70 .seealso: TSGetType() 71 @*/ 72 PetscErrorCode TSSetFromOptions(TS ts) 73 { 74 PetscBool opt,flg; 75 PetscErrorCode ierr; 76 PetscViewer monviewer; 77 char monfilename[PETSC_MAX_PATH_LEN]; 78 SNES snes; 79 TSAdapt adapt; 80 PetscReal time_step; 81 82 PetscFunctionBegin; 83 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 84 ierr = PetscObjectOptionsBegin((PetscObject)ts);CHKERRQ(ierr); 85 /* Handle TS type options */ 86 ierr = TSSetTypeFromOptions(ts);CHKERRQ(ierr); 87 88 /* Handle generic TS options */ 89 ierr = PetscOptionsInt("-ts_max_steps","Maximum number of time steps","TSSetDuration",ts->max_steps,&ts->max_steps,PETSC_NULL);CHKERRQ(ierr); 90 ierr = PetscOptionsReal("-ts_final_time","Time to run to","TSSetDuration",ts->max_time,&ts->max_time,PETSC_NULL);CHKERRQ(ierr); 91 ierr = PetscOptionsReal("-ts_init_time","Initial time","TSSetTime",ts->ptime,&ts->ptime,PETSC_NULL);CHKERRQ(ierr); 92 ierr = PetscOptionsReal("-ts_dt","Initial time step","TSSetTimeStep",ts->time_step,&time_step,&flg);CHKERRQ(ierr); 93 if (flg) { 94 ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 95 } 96 opt = ts->exact_final_time == PETSC_DECIDE ? PETSC_FALSE : (PetscBool)ts->exact_final_time; 97 ierr = PetscOptionsBool("-ts_exact_final_time","Interpolate output to stop exactly at the final time","TSSetExactFinalTime",opt,&opt,&flg);CHKERRQ(ierr); 98 if (flg) {ierr = TSSetExactFinalTime(ts,opt);CHKERRQ(ierr);} 99 ierr = PetscOptionsInt("-ts_max_snes_failures","Maximum number of nonlinear solve failures","TSSetMaxSNESFailures",ts->max_snes_failures,&ts->max_snes_failures,PETSC_NULL);CHKERRQ(ierr); 100 ierr = PetscOptionsInt("-ts_max_reject","Maximum number of step rejections before step fails","TSSetMaxStepRejections",ts->max_reject,&ts->max_reject,PETSC_NULL);CHKERRQ(ierr); 101 ierr = PetscOptionsBool("-ts_error_if_step_fails","Error if no step succeeds","TSSetErrorIfStepFails",ts->errorifstepfailed,&ts->errorifstepfailed,PETSC_NULL);CHKERRQ(ierr); 102 ierr = PetscOptionsReal("-ts_rtol","Relative tolerance for local truncation error","TSSetTolerances",ts->rtol,&ts->rtol,PETSC_NULL);CHKERRQ(ierr); 103 ierr = PetscOptionsReal("-ts_atol","Absolute tolerance for local truncation error","TSSetTolerances",ts->atol,&ts->atol,PETSC_NULL);CHKERRQ(ierr); 104 105 /* Monitor options */ 106 ierr = PetscOptionsString("-ts_monitor","Monitor timestep size","TSMonitorDefault","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 107 if (flg) { 108 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,monfilename,&monviewer);CHKERRQ(ierr); 109 ierr = TSMonitorSet(ts,TSMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 110 } 111 ierr = PetscOptionsString("-ts_monitor_python","Use Python function","TSMonitorSet",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 112 if (flg) {ierr = PetscPythonMonitorSet((PetscObject)ts,monfilename);CHKERRQ(ierr);} 113 114 opt = PETSC_FALSE; 115 ierr = PetscOptionsBool("-ts_monitor_draw","Monitor timestep size graphically","TSMonitorLG",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 116 if (opt) { 117 PetscDrawLG lg; 118 119 ierr = TSMonitorLGCreate(((PetscObject)ts)->comm,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&lg); 120 ierr = TSMonitorSet(ts,TSMonitorLG,lg,TSMonitorLGDestroy);CHKERRQ(ierr); 121 } 122 opt = PETSC_FALSE; 123 ierr = PetscOptionsBool("-ts_monitor_solution","Monitor solution graphically","TSMonitorSolution",opt,&opt,PETSC_NULL);CHKERRQ(ierr); 124 if (opt) { 125 void *ctx; 126 PetscViewer viewer; 127 128 ierr = PetscViewerDrawOpen(((PetscObject)ts)->comm,0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer); 129 ierr = TSMonitorSolutionCreate(ts,viewer,&ctx);CHKERRQ(ierr); 130 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 131 ierr = TSMonitorSet(ts,TSMonitorSolution,ctx,TSMonitorSolutionDestroy);CHKERRQ(ierr); 132 } 133 opt = PETSC_FALSE; 134 ierr = PetscOptionsString("-ts_monitor_solution_binary","Save each solution to a binary file","TSMonitorSolutionBinary",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 135 if (flg) { 136 PetscViewer ctx; 137 if (monfilename[0]) { 138 ierr = PetscViewerBinaryOpen(((PetscObject)ts)->comm,monfilename,FILE_MODE_WRITE,&ctx);CHKERRQ(ierr); 139 } else { 140 ctx = PETSC_VIEWER_BINARY_(((PetscObject)ts)->comm); 141 } 142 ierr = TSMonitorSet(ts,TSMonitorSolutionBinary,ctx,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 143 } 144 opt = PETSC_FALSE; 145 ierr = PetscOptionsString("-ts_monitor_solution_vtk","Save each time step to a binary file, use filename-%%03D.vts","TSMonitorSolutionVTK",0,monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 146 if (flg) { 147 const char *ptr,*ptr2; 148 char *filetemplate; 149 if (!monfilename[0]) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts"); 150 /* Do some cursory validation of the input. */ 151 ierr = PetscStrstr(monfilename,"%",(char**)&ptr);CHKERRQ(ierr); 152 if (!ptr) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"-ts_monitor_solution_vtk requires a file template, e.g. filename-%%03D.vts"); 153 for (ptr++ ; ptr && *ptr; ptr++) { 154 ierr = PetscStrchr("DdiouxX",*ptr,(char**)&ptr2);CHKERRQ(ierr); 155 if (!ptr2 && (*ptr < '0' || '9' < *ptr)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Invalid file template argument to -ts_monitor_solution_vtk, should look like filename-%%03D.vts"); 156 if (ptr2) break; 157 } 158 ierr = PetscStrallocpy(monfilename,&filetemplate);CHKERRQ(ierr); 159 ierr = TSMonitorSet(ts,TSMonitorSolutionVTK,filetemplate,(PetscErrorCode (*)(void**))TSMonitorSolutionVTKDestroy);CHKERRQ(ierr); 160 } 161 162 ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 163 ierr = TSAdaptSetFromOptions(adapt);CHKERRQ(ierr); 164 165 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 166 if (ts->problem_type == TS_LINEAR) {ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr);} 167 168 /* Handle specific TS options */ 169 if (ts->ops->setfromoptions) { 170 ierr = (*ts->ops->setfromoptions)(ts);CHKERRQ(ierr); 171 } 172 173 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 174 ierr = PetscObjectProcessOptionsHandlers((PetscObject)ts);CHKERRQ(ierr); 175 ierr = PetscOptionsEnd();CHKERRQ(ierr); 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #undef __FUNCT__ 181 #define __FUNCT__ "TSComputeRHSJacobian" 182 /*@ 183 TSComputeRHSJacobian - Computes the Jacobian matrix that has been 184 set with TSSetRHSJacobian(). 185 186 Collective on TS and Vec 187 188 Input Parameters: 189 + ts - the TS context 190 . t - current timestep 191 - x - input vector 192 193 Output Parameters: 194 + A - Jacobian matrix 195 . B - optional preconditioning matrix 196 - flag - flag indicating matrix structure 197 198 Notes: 199 Most users should not need to explicitly call this routine, as it 200 is used internally within the nonlinear solvers. 201 202 See KSPSetOperators() for important information about setting the 203 flag parameter. 204 205 Level: developer 206 207 .keywords: SNES, compute, Jacobian, matrix 208 209 .seealso: TSSetRHSJacobian(), KSPSetOperators() 210 @*/ 211 PetscErrorCode TSComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg) 212 { 213 PetscErrorCode ierr; 214 PetscInt Xstate; 215 DM dm; 216 TSDM tsdm; 217 TSRHSJacobian rhsjacobianfunc; 218 void *ctx; 219 TSIJacobian ijacobianfunc; 220 221 PetscFunctionBegin; 222 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 223 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 224 PetscCheckSameComm(ts,1,X,3); 225 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 226 ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr); 227 ierr = DMTSGetRHSJacobian(dm,&rhsjacobianfunc,&ctx);CHKERRQ(ierr); 228 ierr = DMTSGetIJacobian(dm,&ijacobianfunc,PETSC_NULL);CHKERRQ(ierr); 229 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 230 if (ts->rhsjacobian.time == t && (ts->problem_type == TS_LINEAR || (ts->rhsjacobian.X == X && ts->rhsjacobian.Xstate == Xstate))) { 231 *flg = ts->rhsjacobian.mstructure; 232 PetscFunctionReturn(0); 233 } 234 235 if (!rhsjacobianfunc && !ijacobianfunc) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 236 237 if (rhsjacobianfunc) { 238 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 239 *flg = DIFFERENT_NONZERO_PATTERN; 240 PetscStackPush("TS user Jacobian function"); 241 ierr = (*rhsjacobianfunc)(ts,t,X,A,B,flg,ctx);CHKERRQ(ierr); 242 PetscStackPop; 243 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 244 /* make sure user returned a correct Jacobian and preconditioner */ 245 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 246 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 247 } else { 248 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 249 if (*A != *B) {ierr = MatZeroEntries(*B);CHKERRQ(ierr);} 250 *flg = SAME_NONZERO_PATTERN; 251 } 252 ts->rhsjacobian.time = t; 253 ts->rhsjacobian.X = X; 254 ierr = PetscObjectStateQuery((PetscObject)X,&ts->rhsjacobian.Xstate);CHKERRQ(ierr); 255 ts->rhsjacobian.mstructure = *flg; 256 PetscFunctionReturn(0); 257 } 258 259 #undef __FUNCT__ 260 #define __FUNCT__ "TSComputeRHSFunction" 261 /*@ 262 TSComputeRHSFunction - Evaluates the right-hand-side function. 263 264 Collective on TS and Vec 265 266 Input Parameters: 267 + ts - the TS context 268 . t - current time 269 - x - state vector 270 271 Output Parameter: 272 . y - right hand side 273 274 Note: 275 Most users should not need to explicitly call this routine, as it 276 is used internally within the nonlinear solvers. 277 278 Level: developer 279 280 .keywords: TS, compute 281 282 .seealso: TSSetRHSFunction(), TSComputeIFunction() 283 @*/ 284 PetscErrorCode TSComputeRHSFunction(TS ts,PetscReal t,Vec x,Vec y) 285 { 286 PetscErrorCode ierr; 287 288 TSRHSFunction rhsfunction; 289 TSIFunction ifunction; 290 void *ctx; 291 DM dm; 292 293 PetscFunctionBegin; 294 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 295 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 296 PetscValidHeaderSpecific(y,VEC_CLASSID,4); 297 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 298 ierr = DMTSGetRHSFunction(dm,&rhsfunction,&ctx);CHKERRQ(ierr); 299 ierr = DMTSGetIFunction(dm,&ifunction,PETSC_NULL);CHKERRQ(ierr); 300 301 if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 302 303 ierr = PetscLogEventBegin(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 304 if (rhsfunction) { 305 PetscStackPush("TS user right-hand-side function"); 306 ierr = (*rhsfunction)(ts,t,x,y,ctx);CHKERRQ(ierr); 307 PetscStackPop; 308 } else { 309 ierr = VecZeroEntries(y);CHKERRQ(ierr); 310 } 311 312 ierr = PetscLogEventEnd(TS_FunctionEval,ts,x,y,0);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 316 #undef __FUNCT__ 317 #define __FUNCT__ "TSGetRHSVec_Private" 318 static PetscErrorCode TSGetRHSVec_Private(TS ts,Vec *Frhs) 319 { 320 Vec F; 321 PetscErrorCode ierr; 322 323 PetscFunctionBegin; 324 *Frhs = PETSC_NULL; 325 ierr = TSGetIFunction(ts,&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 326 if (!ts->Frhs) { 327 ierr = VecDuplicate(F,&ts->Frhs);CHKERRQ(ierr); 328 } 329 *Frhs = ts->Frhs; 330 PetscFunctionReturn(0); 331 } 332 333 #undef __FUNCT__ 334 #define __FUNCT__ "TSGetRHSMats_Private" 335 static PetscErrorCode TSGetRHSMats_Private(TS ts,Mat *Arhs,Mat *Brhs) 336 { 337 Mat A,B; 338 PetscErrorCode ierr; 339 340 PetscFunctionBegin; 341 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 342 if (Arhs) { 343 if (!ts->Arhs) { 344 ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&ts->Arhs);CHKERRQ(ierr); 345 } 346 *Arhs = ts->Arhs; 347 } 348 if (Brhs) { 349 if (!ts->Brhs) { 350 ierr = MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&ts->Brhs);CHKERRQ(ierr); 351 } 352 *Brhs = ts->Brhs; 353 } 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "TSComputeIFunction" 359 /*@ 360 TSComputeIFunction - Evaluates the DAE residual written in implicit form F(t,X,Xdot)=0 361 362 Collective on TS and Vec 363 364 Input Parameters: 365 + ts - the TS context 366 . t - current time 367 . X - state vector 368 . Xdot - time derivative of state vector 369 - imex - flag indicates if the method is IMEX so that the RHSFunction should be kept separate 370 371 Output Parameter: 372 . Y - right hand side 373 374 Note: 375 Most users should not need to explicitly call this routine, as it 376 is used internally within the nonlinear solvers. 377 378 If the user did did not write their equations in implicit form, this 379 function recasts them in implicit form. 380 381 Level: developer 382 383 .keywords: TS, compute 384 385 .seealso: TSSetIFunction(), TSComputeRHSFunction() 386 @*/ 387 PetscErrorCode TSComputeIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec Y,PetscBool imex) 388 { 389 PetscErrorCode ierr; 390 TSIFunction ifunction; 391 TSRHSFunction rhsfunction; 392 void *ctx; 393 DM dm; 394 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 397 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 398 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 399 PetscValidHeaderSpecific(Y,VEC_CLASSID,5); 400 401 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 402 ierr = DMTSGetIFunction(dm,&ifunction,&ctx);CHKERRQ(ierr); 403 ierr = DMTSGetRHSFunction(dm,&rhsfunction,PETSC_NULL);CHKERRQ(ierr); 404 405 if (!rhsfunction && !ifunction) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSFunction() and / or TSSetIFunction()"); 406 407 ierr = PetscLogEventBegin(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 408 if (ifunction) { 409 PetscStackPush("TS user implicit function"); 410 ierr = (*ifunction)(ts,t,X,Xdot,Y,ctx);CHKERRQ(ierr); 411 PetscStackPop; 412 } 413 if (imex) { 414 if (!ifunction) { 415 ierr = VecCopy(Xdot,Y);CHKERRQ(ierr); 416 } 417 } else if (rhsfunction) { 418 if (ifunction) { 419 Vec Frhs; 420 ierr = TSGetRHSVec_Private(ts,&Frhs);CHKERRQ(ierr); 421 ierr = TSComputeRHSFunction(ts,t,X,Frhs);CHKERRQ(ierr); 422 ierr = VecAXPY(Y,-1,Frhs);CHKERRQ(ierr); 423 } else { 424 ierr = TSComputeRHSFunction(ts,t,X,Y);CHKERRQ(ierr); 425 ierr = VecAYPX(Y,-1,Xdot);CHKERRQ(ierr); 426 } 427 } 428 ierr = PetscLogEventEnd(TS_FunctionEval,ts,X,Xdot,Y);CHKERRQ(ierr); 429 PetscFunctionReturn(0); 430 } 431 432 #undef __FUNCT__ 433 #define __FUNCT__ "TSComputeIJacobian" 434 /*@ 435 TSComputeIJacobian - Evaluates the Jacobian of the DAE 436 437 Collective on TS and Vec 438 439 Input 440 Input Parameters: 441 + ts - the TS context 442 . t - current timestep 443 . X - state vector 444 . Xdot - time derivative of state vector 445 . shift - shift to apply, see note below 446 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 447 448 Output Parameters: 449 + A - Jacobian matrix 450 . B - optional preconditioning matrix 451 - flag - flag indicating matrix structure 452 453 Notes: 454 If F(t,X,Xdot)=0 is the DAE, the required Jacobian is 455 456 dF/dX + shift*dF/dXdot 457 458 Most users should not need to explicitly call this routine, as it 459 is used internally within the nonlinear solvers. 460 461 Level: developer 462 463 .keywords: TS, compute, Jacobian, matrix 464 465 .seealso: TSSetIJacobian() 466 @*/ 467 PetscErrorCode TSComputeIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,PetscBool imex) 468 { 469 PetscInt Xstate, Xdotstate; 470 PetscErrorCode ierr; 471 TSIJacobian ijacobian; 472 TSRHSJacobian rhsjacobian; 473 DM dm; 474 void *ctx; 475 476 PetscFunctionBegin; 477 478 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 479 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 480 PetscValidHeaderSpecific(Xdot,VEC_CLASSID,4); 481 PetscValidPointer(A,6); 482 PetscValidHeaderSpecific(*A,MAT_CLASSID,6); 483 PetscValidPointer(B,7); 484 PetscValidHeaderSpecific(*B,MAT_CLASSID,7); 485 PetscValidPointer(flg,8); 486 487 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 488 ierr = DMTSGetIJacobian(dm,&ijacobian,&ctx);CHKERRQ(ierr); 489 ierr = DMTSGetRHSJacobian(dm,&rhsjacobian,PETSC_NULL);CHKERRQ(ierr); 490 491 ierr = PetscObjectStateQuery((PetscObject)X,&Xstate);CHKERRQ(ierr); 492 ierr = PetscObjectStateQuery((PetscObject)Xdot,&Xdotstate);CHKERRQ(ierr); 493 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))) { 494 *flg = ts->ijacobian.mstructure; 495 ierr = MatScale(*A, shift / ts->ijacobian.shift);CHKERRQ(ierr); 496 PetscFunctionReturn(0); 497 } 498 499 if (!rhsjacobian && !ijacobian) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_USER,"Must call TSSetRHSJacobian() and / or TSSetIJacobian()"); 500 501 *flg = SAME_NONZERO_PATTERN; /* In case we're solving a linear problem in which case it wouldn't get initialized below. */ 502 ierr = PetscLogEventBegin(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 503 if (ijacobian) { 504 *flg = DIFFERENT_NONZERO_PATTERN; 505 PetscStackPush("TS user implicit Jacobian"); 506 ierr = (*ijacobian)(ts,t,X,Xdot,shift,A,B,flg,ctx);CHKERRQ(ierr); 507 PetscStackPop; 508 /* make sure user returned a correct Jacobian and preconditioner */ 509 PetscValidHeaderSpecific(*A,MAT_CLASSID,4); 510 PetscValidHeaderSpecific(*B,MAT_CLASSID,5); 511 } 512 if (imex) { 513 if (!ijacobian) { /* system was written as Xdot = F(t,X) */ 514 ierr = MatZeroEntries(*A);CHKERRQ(ierr); 515 ierr = MatShift(*A,shift);CHKERRQ(ierr); 516 if (*A != *B) { 517 ierr = MatZeroEntries(*B);CHKERRQ(ierr); 518 ierr = MatShift(*B,shift);CHKERRQ(ierr); 519 } 520 *flg = SAME_PRECONDITIONER; 521 } 522 } else { 523 if (!ijacobian) { 524 ierr = TSComputeRHSJacobian(ts,t,X,A,B,flg);CHKERRQ(ierr); 525 ierr = MatScale(*A,-1);CHKERRQ(ierr); 526 ierr = MatShift(*A,shift);CHKERRQ(ierr); 527 if (*A != *B) { 528 ierr = MatScale(*B,-1);CHKERRQ(ierr); 529 ierr = MatShift(*B,shift);CHKERRQ(ierr); 530 } 531 } else if (rhsjacobian) { 532 Mat Arhs,Brhs; 533 MatStructure axpy,flg2 = DIFFERENT_NONZERO_PATTERN; 534 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 535 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 536 axpy = (*flg == flg2) ? SAME_NONZERO_PATTERN : DIFFERENT_NONZERO_PATTERN; 537 ierr = MatAXPY(*A,-1,Arhs,axpy);CHKERRQ(ierr); 538 if (*A != *B) { 539 ierr = MatAXPY(*B,-1,Brhs,axpy);CHKERRQ(ierr); 540 } 541 *flg = PetscMin(*flg,flg2); 542 } 543 } 544 545 ts->ijacobian.time = t; 546 ts->ijacobian.X = X; 547 ts->ijacobian.Xdot = Xdot; 548 ierr = PetscObjectStateQuery((PetscObject)X,&ts->ijacobian.Xstate);CHKERRQ(ierr); 549 ierr = PetscObjectStateQuery((PetscObject)Xdot,&ts->ijacobian.Xdotstate);CHKERRQ(ierr); 550 ts->ijacobian.shift = shift; 551 ts->ijacobian.imex = imex; 552 ts->ijacobian.mstructure = *flg; 553 ierr = PetscLogEventEnd(TS_JacobianEval,ts,X,*A,*B);CHKERRQ(ierr); 554 PetscFunctionReturn(0); 555 } 556 557 #undef __FUNCT__ 558 #define __FUNCT__ "TSSetRHSFunction" 559 /*@C 560 TSSetRHSFunction - Sets the routine for evaluating the function, 561 F(t,u), where U_t = F(t,u). 562 563 Logically Collective on TS 564 565 Input Parameters: 566 + ts - the TS context obtained from TSCreate() 567 . r - vector to put the computed right hand side (or PETSC_NULL to have it created) 568 . f - routine for evaluating the right-hand-side function 569 - ctx - [optional] user-defined context for private data for the 570 function evaluation routine (may be PETSC_NULL) 571 572 Calling sequence of func: 573 $ func (TS ts,PetscReal t,Vec u,Vec F,void *ctx); 574 575 + t - current timestep 576 . u - input vector 577 . F - function vector 578 - ctx - [optional] user-defined function context 579 580 Level: beginner 581 582 .keywords: TS, timestep, set, right-hand-side, function 583 584 .seealso: TSSetRHSJacobian(), TSSetIJacobian() 585 @*/ 586 PetscErrorCode TSSetRHSFunction(TS ts,Vec r,PetscErrorCode (*f)(TS,PetscReal,Vec,Vec,void*),void *ctx) 587 { 588 PetscErrorCode ierr; 589 SNES snes; 590 Vec ralloc = PETSC_NULL; 591 DM dm; 592 593 PetscFunctionBegin; 594 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 595 if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 596 597 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 598 ierr = DMTSSetRHSFunction(dm,f,ctx);CHKERRQ(ierr); 599 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 600 if (!r && !ts->dm && ts->vec_sol) { 601 ierr = VecDuplicate(ts->vec_sol,&ralloc);CHKERRQ(ierr); 602 r = ralloc; 603 } 604 ierr = SNESSetFunction(snes,r,SNESTSFormFunction,ts);CHKERRQ(ierr); 605 ierr = VecDestroy(&ralloc);CHKERRQ(ierr); 606 PetscFunctionReturn(0); 607 } 608 609 #undef __FUNCT__ 610 #define __FUNCT__ "TSSetRHSJacobian" 611 /*@C 612 TSSetRHSJacobian - Sets the function to compute the Jacobian of F, 613 where U_t = F(U,t), as well as the location to store the matrix. 614 615 Logically Collective on TS 616 617 Input Parameters: 618 + ts - the TS context obtained from TSCreate() 619 . A - Jacobian matrix 620 . B - preconditioner matrix (usually same as A) 621 . f - the Jacobian evaluation routine 622 - ctx - [optional] user-defined context for private data for the 623 Jacobian evaluation routine (may be PETSC_NULL) 624 625 Calling sequence of func: 626 $ func (TS ts,PetscReal t,Vec u,Mat *A,Mat *B,MatStructure *flag,void *ctx); 627 628 + t - current timestep 629 . u - input vector 630 . A - matrix A, where U_t = A(t)u 631 . B - preconditioner matrix, usually the same as A 632 . flag - flag indicating information about the preconditioner matrix 633 structure (same as flag in KSPSetOperators()) 634 - ctx - [optional] user-defined context for matrix evaluation routine 635 636 Notes: 637 See KSPSetOperators() for important information about setting the flag 638 output parameter in the routine func(). Be sure to read this information! 639 640 The routine func() takes Mat * as the matrix arguments rather than Mat. 641 This allows the matrix evaluation routine to replace A and/or B with a 642 completely new matrix structure (not just different matrix elements) 643 when appropriate, for instance, if the nonzero structure is changing 644 throughout the global iterations. 645 646 Level: beginner 647 648 .keywords: TS, timestep, set, right-hand-side, Jacobian 649 650 .seealso: SNESDefaultComputeJacobianColor(), TSSetRHSFunction() 651 652 @*/ 653 PetscErrorCode TSSetRHSJacobian(TS ts,Mat A,Mat B,TSRHSJacobian f,void *ctx) 654 { 655 PetscErrorCode ierr; 656 SNES snes; 657 DM dm; 658 TSIJacobian ijacobian; 659 660 PetscFunctionBegin; 661 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 662 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 663 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 664 if (A) PetscCheckSameComm(ts,1,A,2); 665 if (B) PetscCheckSameComm(ts,1,B,3); 666 667 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 668 ierr = DMTSSetRHSJacobian(dm,f,ctx);CHKERRQ(ierr); 669 ierr = DMTSGetIJacobian(dm,&ijacobian,PETSC_NULL);CHKERRQ(ierr); 670 671 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 672 if (!ijacobian) { 673 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 674 } 675 if (A) { 676 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 677 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 678 ts->Arhs = A; 679 } 680 if (B) { 681 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 682 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 683 ts->Brhs = B; 684 } 685 PetscFunctionReturn(0); 686 } 687 688 689 #undef __FUNCT__ 690 #define __FUNCT__ "TSSetIFunction" 691 /*@C 692 TSSetIFunction - Set the function to compute F(t,U,U_t) where F = 0 is the DAE to be solved. 693 694 Logically Collective on TS 695 696 Input Parameters: 697 + ts - the TS context obtained from TSCreate() 698 . r - vector to hold the residual (or PETSC_NULL to have it created internally) 699 . f - the function evaluation routine 700 - ctx - user-defined context for private data for the function evaluation routine (may be PETSC_NULL) 701 702 Calling sequence of f: 703 $ f(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 704 705 + t - time at step/stage being solved 706 . u - state vector 707 . u_t - time derivative of state vector 708 . F - function vector 709 - ctx - [optional] user-defined context for matrix evaluation routine 710 711 Important: 712 The user MUST call either this routine, TSSetRHSFunction(). This routine must be used when not solving an ODE, for example a DAE. 713 714 Level: beginner 715 716 .keywords: TS, timestep, set, DAE, Jacobian 717 718 .seealso: TSSetRHSJacobian(), TSSetRHSFunction(), TSSetIJacobian() 719 @*/ 720 PetscErrorCode TSSetIFunction(TS ts,Vec res,TSIFunction f,void *ctx) 721 { 722 PetscErrorCode ierr; 723 SNES snes; 724 Vec resalloc = PETSC_NULL; 725 DM dm; 726 727 PetscFunctionBegin; 728 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 729 if (res) PetscValidHeaderSpecific(res,VEC_CLASSID,2); 730 731 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 732 ierr = DMTSSetIFunction(dm,f,ctx);CHKERRQ(ierr); 733 734 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 735 if (!res && !ts->dm && ts->vec_sol) { 736 ierr = VecDuplicate(ts->vec_sol,&resalloc);CHKERRQ(ierr); 737 res = resalloc; 738 } 739 ierr = SNESSetFunction(snes,res,SNESTSFormFunction,ts);CHKERRQ(ierr); 740 ierr = VecDestroy(&resalloc);CHKERRQ(ierr); 741 742 PetscFunctionReturn(0); 743 } 744 745 #undef __FUNCT__ 746 #define __FUNCT__ "TSGetIFunction" 747 /*@C 748 TSGetIFunction - Returns the vector where the implicit residual is stored and the function/contex to compute it. 749 750 Not Collective 751 752 Input Parameter: 753 . ts - the TS context 754 755 Output Parameter: 756 + r - vector to hold residual (or PETSC_NULL) 757 . func - the function to compute residual (or PETSC_NULL) 758 - ctx - the function context (or PETSC_NULL) 759 760 Level: advanced 761 762 .keywords: TS, nonlinear, get, function 763 764 .seealso: TSSetIFunction(), SNESGetFunction() 765 @*/ 766 PetscErrorCode TSGetIFunction(TS ts,Vec *r,TSIFunction *func,void **ctx) 767 { 768 PetscErrorCode ierr; 769 SNES snes; 770 DM dm; 771 772 PetscFunctionBegin; 773 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 774 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 775 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 776 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 777 ierr = DMTSGetIFunction(dm,func,ctx);CHKERRQ(ierr); 778 PetscFunctionReturn(0); 779 } 780 781 #undef __FUNCT__ 782 #define __FUNCT__ "TSGetRHSFunction" 783 /*@C 784 TSGetRHSFunction - Returns the vector where the right hand side is stored and the function/context to compute it. 785 786 Not Collective 787 788 Input Parameter: 789 . ts - the TS context 790 791 Output Parameter: 792 + r - vector to hold computed right hand side (or PETSC_NULL) 793 . func - the function to compute right hand side (or PETSC_NULL) 794 - ctx - the function context (or PETSC_NULL) 795 796 Level: advanced 797 798 .keywords: TS, nonlinear, get, function 799 800 .seealso: TSSetRhsfunction(), SNESGetFunction() 801 @*/ 802 PetscErrorCode TSGetRHSFunction(TS ts,Vec *r,TSRHSFunction *func,void **ctx) 803 { 804 PetscErrorCode ierr; 805 SNES snes; 806 DM dm; 807 808 PetscFunctionBegin; 809 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 810 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 811 ierr = SNESGetFunction(snes,r,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 812 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 813 ierr = DMTSGetRHSFunction(dm,func,ctx);CHKERRQ(ierr); 814 PetscFunctionReturn(0); 815 } 816 817 #undef __FUNCT__ 818 #define __FUNCT__ "TSSetIJacobian" 819 /*@C 820 TSSetIJacobian - Set the function to compute the matrix dF/dU + a*dF/dU_t where F(t,U,U_t) is the function 821 you provided with TSSetIFunction(). 822 823 Logically Collective on TS 824 825 Input Parameters: 826 + ts - the TS context obtained from TSCreate() 827 . A - Jacobian matrix 828 . B - preconditioning matrix for A (may be same as A) 829 . f - the Jacobian evaluation routine 830 - ctx - user-defined context for private data for the Jacobian evaluation routine (may be PETSC_NULL) 831 832 Calling sequence of f: 833 $ f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat *A,Mat *B,MatStructure *flag,void *ctx); 834 835 + t - time at step/stage being solved 836 . U - state vector 837 . U_t - time derivative of state vector 838 . a - shift 839 . A - Jacobian of G(U) = F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 840 . B - preconditioning matrix for A, may be same as A 841 . flag - flag indicating information about the preconditioner matrix 842 structure (same as flag in KSPSetOperators()) 843 - ctx - [optional] user-defined context for matrix evaluation routine 844 845 Notes: 846 The matrices A and B are exactly the matrices that are used by SNES for the nonlinear solve. 847 848 The matrix dF/dU + a*dF/dU_t you provide turns out to be 849 the Jacobian of G(U) = F(t,U,W+a*U) where F(t,U,U_t) = 0 is the DAE to be solved. 850 The time integrator internally approximates U_t by W+a*U where the positive "shift" 851 a and vector W depend on the integration method, step size, and past states. For example with 852 the backward Euler method a = 1/dt and W = -a*U(previous timestep) so 853 W + a*U = a*(U - U(previous timestep)) = (U - U(previous timestep))/dt 854 855 Level: beginner 856 857 .keywords: TS, timestep, DAE, Jacobian 858 859 .seealso: TSSetIFunction(), TSSetRHSJacobian(), SNESDefaultComputeJacobianColor(), SNESDefaultComputeJacobian() 860 861 @*/ 862 PetscErrorCode TSSetIJacobian(TS ts,Mat A,Mat B,TSIJacobian f,void *ctx) 863 { 864 PetscErrorCode ierr; 865 SNES snes; 866 DM dm; 867 868 PetscFunctionBegin; 869 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 870 if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,2); 871 if (B) PetscValidHeaderSpecific(B,MAT_CLASSID,3); 872 if (A) PetscCheckSameComm(ts,1,A,2); 873 if (B) PetscCheckSameComm(ts,1,B,3); 874 875 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 876 ierr = DMTSSetIJacobian(dm,f,ctx);CHKERRQ(ierr); 877 878 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 879 ierr = SNESSetJacobian(snes,A,B,SNESTSFormJacobian,ts);CHKERRQ(ierr); 880 PetscFunctionReturn(0); 881 } 882 883 #undef __FUNCT__ 884 #define __FUNCT__ "TSView" 885 /*@C 886 TSView - Prints the TS data structure. 887 888 Collective on TS 889 890 Input Parameters: 891 + ts - the TS context obtained from TSCreate() 892 - viewer - visualization context 893 894 Options Database Key: 895 . -ts_view - calls TSView() at end of TSStep() 896 897 Notes: 898 The available visualization contexts include 899 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 900 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 901 output where only the first processor opens 902 the file. All other processors send their 903 data to the first processor to print. 904 905 The user can open an alternative visualization context with 906 PetscViewerASCIIOpen() - output to a specified file. 907 908 Level: beginner 909 910 .keywords: TS, timestep, view 911 912 .seealso: PetscViewerASCIIOpen() 913 @*/ 914 PetscErrorCode TSView(TS ts,PetscViewer viewer) 915 { 916 PetscErrorCode ierr; 917 TSType type; 918 PetscBool iascii,isstring,isundials; 919 920 PetscFunctionBegin; 921 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 922 if (!viewer) { 923 ierr = PetscViewerASCIIGetStdout(((PetscObject)ts)->comm,&viewer);CHKERRQ(ierr); 924 } 925 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 926 PetscCheckSameComm(ts,1,viewer,2); 927 928 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 929 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 930 if (iascii) { 931 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)ts,viewer,"TS Object");CHKERRQ(ierr); 932 ierr = PetscViewerASCIIPrintf(viewer," maximum steps=%D\n",ts->max_steps);CHKERRQ(ierr); 933 ierr = PetscViewerASCIIPrintf(viewer," maximum time=%G\n",ts->max_time);CHKERRQ(ierr); 934 if (ts->problem_type == TS_NONLINEAR) { 935 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solver iterations=%D\n",ts->snes_its);CHKERRQ(ierr); 936 ierr = PetscViewerASCIIPrintf(viewer," total number of nonlinear solve failures=%D\n",ts->num_snes_failures);CHKERRQ(ierr); 937 } 938 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",ts->ksp_its);CHKERRQ(ierr); 939 ierr = PetscViewerASCIIPrintf(viewer," total number of rejected steps=%D\n",ts->reject);CHKERRQ(ierr); 940 if (ts->ops->view) { 941 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 942 ierr = (*ts->ops->view)(ts,viewer);CHKERRQ(ierr); 943 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 944 } 945 } else if (isstring) { 946 ierr = TSGetType(ts,&type);CHKERRQ(ierr); 947 ierr = PetscViewerStringSPrintf(viewer," %-7.7s",type);CHKERRQ(ierr); 948 } 949 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 950 ierr = PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&isundials);CHKERRQ(ierr); 951 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 952 PetscFunctionReturn(0); 953 } 954 955 956 #undef __FUNCT__ 957 #define __FUNCT__ "TSSetApplicationContext" 958 /*@ 959 TSSetApplicationContext - Sets an optional user-defined context for 960 the timesteppers. 961 962 Logically Collective on TS 963 964 Input Parameters: 965 + ts - the TS context obtained from TSCreate() 966 - usrP - optional user context 967 968 Level: intermediate 969 970 .keywords: TS, timestep, set, application, context 971 972 .seealso: TSGetApplicationContext() 973 @*/ 974 PetscErrorCode TSSetApplicationContext(TS ts,void *usrP) 975 { 976 PetscFunctionBegin; 977 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 978 ts->user = usrP; 979 PetscFunctionReturn(0); 980 } 981 982 #undef __FUNCT__ 983 #define __FUNCT__ "TSGetApplicationContext" 984 /*@ 985 TSGetApplicationContext - Gets the user-defined context for the 986 timestepper. 987 988 Not Collective 989 990 Input Parameter: 991 . ts - the TS context obtained from TSCreate() 992 993 Output Parameter: 994 . usrP - user context 995 996 Level: intermediate 997 998 .keywords: TS, timestep, get, application, context 999 1000 .seealso: TSSetApplicationContext() 1001 @*/ 1002 PetscErrorCode TSGetApplicationContext(TS ts,void *usrP) 1003 { 1004 PetscFunctionBegin; 1005 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1006 *(void**)usrP = ts->user; 1007 PetscFunctionReturn(0); 1008 } 1009 1010 #undef __FUNCT__ 1011 #define __FUNCT__ "TSGetTimeStepNumber" 1012 /*@ 1013 TSGetTimeStepNumber - Gets the number of time steps completed. 1014 1015 Not Collective 1016 1017 Input Parameter: 1018 . ts - the TS context obtained from TSCreate() 1019 1020 Output Parameter: 1021 . iter - number of steps completed so far 1022 1023 Level: intermediate 1024 1025 .keywords: TS, timestep, get, iteration, number 1026 .seealso: TSGetTime(), TSGetTimeStep(), TSSetPreStep(), TSSetPreStage(), TSSetPostStep() 1027 @*/ 1028 PetscErrorCode TSGetTimeStepNumber(TS ts,PetscInt* iter) 1029 { 1030 PetscFunctionBegin; 1031 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1032 PetscValidIntPointer(iter,2); 1033 *iter = ts->steps; 1034 PetscFunctionReturn(0); 1035 } 1036 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "TSSetInitialTimeStep" 1039 /*@ 1040 TSSetInitialTimeStep - Sets the initial timestep to be used, 1041 as well as the initial time. 1042 1043 Logically Collective on TS 1044 1045 Input Parameters: 1046 + ts - the TS context obtained from TSCreate() 1047 . initial_time - the initial time 1048 - time_step - the size of the timestep 1049 1050 Level: intermediate 1051 1052 .seealso: TSSetTimeStep(), TSGetTimeStep() 1053 1054 .keywords: TS, set, initial, timestep 1055 @*/ 1056 PetscErrorCode TSSetInitialTimeStep(TS ts,PetscReal initial_time,PetscReal time_step) 1057 { 1058 PetscErrorCode ierr; 1059 1060 PetscFunctionBegin; 1061 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1062 ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr); 1063 ierr = TSSetTime(ts,initial_time);CHKERRQ(ierr); 1064 PetscFunctionReturn(0); 1065 } 1066 1067 #undef __FUNCT__ 1068 #define __FUNCT__ "TSSetTimeStep" 1069 /*@ 1070 TSSetTimeStep - Allows one to reset the timestep at any time, 1071 useful for simple pseudo-timestepping codes. 1072 1073 Logically Collective on TS 1074 1075 Input Parameters: 1076 + ts - the TS context obtained from TSCreate() 1077 - time_step - the size of the timestep 1078 1079 Level: intermediate 1080 1081 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1082 1083 .keywords: TS, set, timestep 1084 @*/ 1085 PetscErrorCode TSSetTimeStep(TS ts,PetscReal time_step) 1086 { 1087 PetscFunctionBegin; 1088 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1089 PetscValidLogicalCollectiveReal(ts,time_step,2); 1090 ts->time_step = time_step; 1091 ts->time_step_orig = time_step; 1092 PetscFunctionReturn(0); 1093 } 1094 1095 #undef __FUNCT__ 1096 #define __FUNCT__ "TSSetExactFinalTime" 1097 /*@ 1098 TSSetExactFinalTime - Determines whether to interpolate solution to the 1099 exact final time requested by the user or just returns it at the final time 1100 it computed. 1101 1102 Logically Collective on TS 1103 1104 Input Parameter: 1105 + ts - the time-step context 1106 - ft - PETSC_TRUE if interpolates, else PETSC_FALSE 1107 1108 Level: beginner 1109 1110 .seealso: TSSetDuration() 1111 @*/ 1112 PetscErrorCode TSSetExactFinalTime(TS ts,PetscBool flg) 1113 { 1114 1115 PetscFunctionBegin; 1116 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1117 PetscValidLogicalCollectiveBool(ts,flg,2); 1118 ts->exact_final_time = flg; 1119 PetscFunctionReturn(0); 1120 } 1121 1122 #undef __FUNCT__ 1123 #define __FUNCT__ "TSGetTimeStep" 1124 /*@ 1125 TSGetTimeStep - Gets the current timestep size. 1126 1127 Not Collective 1128 1129 Input Parameter: 1130 . ts - the TS context obtained from TSCreate() 1131 1132 Output Parameter: 1133 . dt - the current timestep size 1134 1135 Level: intermediate 1136 1137 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 1138 1139 .keywords: TS, get, timestep 1140 @*/ 1141 PetscErrorCode TSGetTimeStep(TS ts,PetscReal* dt) 1142 { 1143 PetscFunctionBegin; 1144 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1145 PetscValidRealPointer(dt,2); 1146 *dt = ts->time_step; 1147 PetscFunctionReturn(0); 1148 } 1149 1150 #undef __FUNCT__ 1151 #define __FUNCT__ "TSGetSolution" 1152 /*@ 1153 TSGetSolution - Returns the solution at the present timestep. It 1154 is valid to call this routine inside the function that you are evaluating 1155 in order to move to the new timestep. This vector not changed until 1156 the solution at the next timestep has been calculated. 1157 1158 Not Collective, but Vec returned is parallel if TS is parallel 1159 1160 Input Parameter: 1161 . ts - the TS context obtained from TSCreate() 1162 1163 Output Parameter: 1164 . v - the vector containing the solution 1165 1166 Level: intermediate 1167 1168 .seealso: TSGetTimeStep() 1169 1170 .keywords: TS, timestep, get, solution 1171 @*/ 1172 PetscErrorCode TSGetSolution(TS ts,Vec *v) 1173 { 1174 PetscFunctionBegin; 1175 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1176 PetscValidPointer(v,2); 1177 *v = ts->vec_sol; 1178 PetscFunctionReturn(0); 1179 } 1180 1181 /* ----- Routines to initialize and destroy a timestepper ---- */ 1182 #undef __FUNCT__ 1183 #define __FUNCT__ "TSSetProblemType" 1184 /*@ 1185 TSSetProblemType - Sets the type of problem to be solved. 1186 1187 Not collective 1188 1189 Input Parameters: 1190 + ts - The TS 1191 - type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1192 .vb 1193 U_t = A U 1194 U_t = A(t) U 1195 U_t = F(t,U) 1196 .ve 1197 1198 Level: beginner 1199 1200 .keywords: TS, problem type 1201 .seealso: TSSetUp(), TSProblemType, TS 1202 @*/ 1203 PetscErrorCode TSSetProblemType(TS ts, TSProblemType type) 1204 { 1205 PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1209 ts->problem_type = type; 1210 if (type == TS_LINEAR) { 1211 SNES snes; 1212 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1213 ierr = SNESSetType(snes,SNESKSPONLY);CHKERRQ(ierr); 1214 } 1215 PetscFunctionReturn(0); 1216 } 1217 1218 #undef __FUNCT__ 1219 #define __FUNCT__ "TSGetProblemType" 1220 /*@C 1221 TSGetProblemType - Gets the type of problem to be solved. 1222 1223 Not collective 1224 1225 Input Parameter: 1226 . ts - The TS 1227 1228 Output Parameter: 1229 . type - One of TS_LINEAR, TS_NONLINEAR where these types refer to problems of the forms 1230 .vb 1231 M U_t = A U 1232 M(t) U_t = A(t) U 1233 U_t = F(t,U) 1234 .ve 1235 1236 Level: beginner 1237 1238 .keywords: TS, problem type 1239 .seealso: TSSetUp(), TSProblemType, TS 1240 @*/ 1241 PetscErrorCode TSGetProblemType(TS ts, TSProblemType *type) 1242 { 1243 PetscFunctionBegin; 1244 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1245 PetscValidIntPointer(type,2); 1246 *type = ts->problem_type; 1247 PetscFunctionReturn(0); 1248 } 1249 1250 #undef __FUNCT__ 1251 #define __FUNCT__ "TSSetUp" 1252 /*@ 1253 TSSetUp - Sets up the internal data structures for the later use 1254 of a timestepper. 1255 1256 Collective on TS 1257 1258 Input Parameter: 1259 . ts - the TS context obtained from TSCreate() 1260 1261 Notes: 1262 For basic use of the TS solvers the user need not explicitly call 1263 TSSetUp(), since these actions will automatically occur during 1264 the call to TSStep(). However, if one wishes to control this 1265 phase separately, TSSetUp() should be called after TSCreate() 1266 and optional routines of the form TSSetXXX(), but before TSStep(). 1267 1268 Level: advanced 1269 1270 .keywords: TS, timestep, setup 1271 1272 .seealso: TSCreate(), TSStep(), TSDestroy() 1273 @*/ 1274 PetscErrorCode TSSetUp(TS ts) 1275 { 1276 PetscErrorCode ierr; 1277 DM dm; 1278 PetscErrorCode (*func)(SNES,Vec,Vec,void*); 1279 PetscErrorCode (*jac)(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 1280 TSIJacobian ijac; 1281 TSRHSJacobian rhsjac; 1282 1283 PetscFunctionBegin; 1284 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1285 if (ts->setupcalled) PetscFunctionReturn(0); 1286 1287 if (!((PetscObject)ts)->type_name) { 1288 ierr = TSSetType(ts,TSEULER);CHKERRQ(ierr); 1289 } 1290 if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE; 1291 1292 if (!ts->vec_sol) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetSolution() first"); 1293 1294 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 1295 1296 if (ts->ops->setup) { 1297 ierr = (*ts->ops->setup)(ts);CHKERRQ(ierr); 1298 } 1299 1300 /* in the case where we've set a DMTSFunction or what have you, we need the default SNESFunction 1301 to be set right but can't do it elsewhere due to the overreliance on ctx=ts. 1302 */ 1303 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1304 ierr = DMSNESGetFunction(dm,&func,PETSC_NULL);CHKERRQ(ierr); 1305 if (!func) { 1306 ierr =DMSNESSetFunction(dm,SNESTSFormFunction,ts);CHKERRQ(ierr); 1307 } 1308 /* if the SNES doesn't have a jacobian set and the TS has an ijacobian or rhsjacobian set, set the SNES to use it. 1309 Otherwise, the SNES will use coloring internally to form the Jacobian. 1310 */ 1311 ierr = DMSNESGetJacobian(dm,&jac,PETSC_NULL);CHKERRQ(ierr); 1312 ierr = DMTSGetIJacobian(dm,&ijac,PETSC_NULL);CHKERRQ(ierr); 1313 ierr = DMTSGetRHSJacobian(dm,&rhsjac,PETSC_NULL);CHKERRQ(ierr); 1314 if (!jac && (ijac || rhsjac)) { 1315 ierr = DMSNESSetJacobian(dm,SNESTSFormJacobian,ts);CHKERRQ(ierr); 1316 } 1317 ts->setupcalled = PETSC_TRUE; 1318 PetscFunctionReturn(0); 1319 } 1320 1321 #undef __FUNCT__ 1322 #define __FUNCT__ "TSReset" 1323 /*@ 1324 TSReset - Resets a TS context and removes any allocated Vecs and Mats. 1325 1326 Collective on TS 1327 1328 Input Parameter: 1329 . ts - the TS context obtained from TSCreate() 1330 1331 Level: beginner 1332 1333 .keywords: TS, timestep, reset 1334 1335 .seealso: TSCreate(), TSSetup(), TSDestroy() 1336 @*/ 1337 PetscErrorCode TSReset(TS ts) 1338 { 1339 PetscErrorCode ierr; 1340 1341 PetscFunctionBegin; 1342 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1343 if (ts->ops->reset) { 1344 ierr = (*ts->ops->reset)(ts);CHKERRQ(ierr); 1345 } 1346 if (ts->snes) {ierr = SNESReset(ts->snes);CHKERRQ(ierr);} 1347 ierr = MatDestroy(&ts->Arhs);CHKERRQ(ierr); 1348 ierr = MatDestroy(&ts->Brhs);CHKERRQ(ierr); 1349 ierr = VecDestroy(&ts->Frhs);CHKERRQ(ierr); 1350 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1351 ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr); 1352 ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr); 1353 ierr = VecDestroyVecs(ts->nwork,&ts->work);CHKERRQ(ierr); 1354 ts->setupcalled = PETSC_FALSE; 1355 PetscFunctionReturn(0); 1356 } 1357 1358 #undef __FUNCT__ 1359 #define __FUNCT__ "TSDestroy" 1360 /*@ 1361 TSDestroy - Destroys the timestepper context that was created 1362 with TSCreate(). 1363 1364 Collective on TS 1365 1366 Input Parameter: 1367 . ts - the TS context obtained from TSCreate() 1368 1369 Level: beginner 1370 1371 .keywords: TS, timestepper, destroy 1372 1373 .seealso: TSCreate(), TSSetUp(), TSSolve() 1374 @*/ 1375 PetscErrorCode TSDestroy(TS *ts) 1376 { 1377 PetscErrorCode ierr; 1378 1379 PetscFunctionBegin; 1380 if (!*ts) PetscFunctionReturn(0); 1381 PetscValidHeaderSpecific((*ts),TS_CLASSID,1); 1382 if (--((PetscObject)(*ts))->refct > 0) {*ts = 0; PetscFunctionReturn(0);} 1383 1384 ierr = TSReset((*ts));CHKERRQ(ierr); 1385 1386 /* if memory was published with AMS then destroy it */ 1387 ierr = PetscObjectDepublish((*ts));CHKERRQ(ierr); 1388 if ((*ts)->ops->destroy) {ierr = (*(*ts)->ops->destroy)((*ts));CHKERRQ(ierr);} 1389 1390 ierr = TSAdaptDestroy(&(*ts)->adapt);CHKERRQ(ierr); 1391 ierr = SNESDestroy(&(*ts)->snes);CHKERRQ(ierr); 1392 ierr = DMDestroy(&(*ts)->dm);CHKERRQ(ierr); 1393 ierr = TSMonitorCancel((*ts));CHKERRQ(ierr); 1394 1395 ierr = PetscHeaderDestroy(ts);CHKERRQ(ierr); 1396 PetscFunctionReturn(0); 1397 } 1398 1399 #undef __FUNCT__ 1400 #define __FUNCT__ "TSGetSNES" 1401 /*@ 1402 TSGetSNES - Returns the SNES (nonlinear solver) associated with 1403 a TS (timestepper) context. Valid only for nonlinear problems. 1404 1405 Not Collective, but SNES is parallel if TS is parallel 1406 1407 Input Parameter: 1408 . ts - the TS context obtained from TSCreate() 1409 1410 Output Parameter: 1411 . snes - the nonlinear solver context 1412 1413 Notes: 1414 The user can then directly manipulate the SNES context to set various 1415 options, etc. Likewise, the user can then extract and manipulate the 1416 KSP, KSP, and PC contexts as well. 1417 1418 TSGetSNES() does not work for integrators that do not use SNES; in 1419 this case TSGetSNES() returns PETSC_NULL in snes. 1420 1421 Level: beginner 1422 1423 .keywords: timestep, get, SNES 1424 @*/ 1425 PetscErrorCode TSGetSNES(TS ts,SNES *snes) 1426 { 1427 PetscErrorCode ierr; 1428 1429 PetscFunctionBegin; 1430 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1431 PetscValidPointer(snes,2); 1432 if (!ts->snes) { 1433 ierr = SNESCreate(((PetscObject)ts)->comm,&ts->snes);CHKERRQ(ierr); 1434 ierr = SNESSetFunction(ts->snes,PETSC_NULL,SNESTSFormFunction,ts);CHKERRQ(ierr); 1435 ierr = PetscLogObjectParent(ts,ts->snes);CHKERRQ(ierr); 1436 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->snes,(PetscObject)ts,1);CHKERRQ(ierr); 1437 if (ts->dm) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);} 1438 if (ts->problem_type == TS_LINEAR) { 1439 ierr = SNESSetType(ts->snes,SNESKSPONLY);CHKERRQ(ierr); 1440 } 1441 } 1442 *snes = ts->snes; 1443 PetscFunctionReturn(0); 1444 } 1445 1446 #undef __FUNCT__ 1447 #define __FUNCT__ "TSGetKSP" 1448 /*@ 1449 TSGetKSP - Returns the KSP (linear solver) associated with 1450 a TS (timestepper) context. 1451 1452 Not Collective, but KSP is parallel if TS is parallel 1453 1454 Input Parameter: 1455 . ts - the TS context obtained from TSCreate() 1456 1457 Output Parameter: 1458 . ksp - the nonlinear solver context 1459 1460 Notes: 1461 The user can then directly manipulate the KSP context to set various 1462 options, etc. Likewise, the user can then extract and manipulate the 1463 KSP and PC contexts as well. 1464 1465 TSGetKSP() does not work for integrators that do not use KSP; 1466 in this case TSGetKSP() returns PETSC_NULL in ksp. 1467 1468 Level: beginner 1469 1470 .keywords: timestep, get, KSP 1471 @*/ 1472 PetscErrorCode TSGetKSP(TS ts,KSP *ksp) 1473 { 1474 PetscErrorCode ierr; 1475 SNES snes; 1476 1477 PetscFunctionBegin; 1478 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1479 PetscValidPointer(ksp,2); 1480 if (!((PetscObject)ts)->type_name) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"KSP is not created yet. Call TSSetType() first"); 1481 if (ts->problem_type != TS_LINEAR) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Linear only; use TSGetSNES()"); 1482 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 1483 ierr = SNESGetKSP(snes,ksp);CHKERRQ(ierr); 1484 PetscFunctionReturn(0); 1485 } 1486 1487 /* ----------- Routines to set solver parameters ---------- */ 1488 1489 #undef __FUNCT__ 1490 #define __FUNCT__ "TSGetDuration" 1491 /*@ 1492 TSGetDuration - Gets the maximum number of timesteps to use and 1493 maximum time for iteration. 1494 1495 Not Collective 1496 1497 Input Parameters: 1498 + ts - the TS context obtained from TSCreate() 1499 . maxsteps - maximum number of iterations to use, or PETSC_NULL 1500 - maxtime - final time to iterate to, or PETSC_NULL 1501 1502 Level: intermediate 1503 1504 .keywords: TS, timestep, get, maximum, iterations, time 1505 @*/ 1506 PetscErrorCode TSGetDuration(TS ts, PetscInt *maxsteps, PetscReal *maxtime) 1507 { 1508 PetscFunctionBegin; 1509 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1510 if (maxsteps) { 1511 PetscValidIntPointer(maxsteps,2); 1512 *maxsteps = ts->max_steps; 1513 } 1514 if (maxtime) { 1515 PetscValidScalarPointer(maxtime,3); 1516 *maxtime = ts->max_time; 1517 } 1518 PetscFunctionReturn(0); 1519 } 1520 1521 #undef __FUNCT__ 1522 #define __FUNCT__ "TSSetDuration" 1523 /*@ 1524 TSSetDuration - Sets the maximum number of timesteps to use and 1525 maximum time for iteration. 1526 1527 Logically Collective on TS 1528 1529 Input Parameters: 1530 + ts - the TS context obtained from TSCreate() 1531 . maxsteps - maximum number of iterations to use 1532 - maxtime - final time to iterate to 1533 1534 Options Database Keys: 1535 . -ts_max_steps <maxsteps> - Sets maxsteps 1536 . -ts_final_time <maxtime> - Sets maxtime 1537 1538 Notes: 1539 The default maximum number of iterations is 5000. Default time is 5.0 1540 1541 Level: intermediate 1542 1543 .keywords: TS, timestep, set, maximum, iterations 1544 1545 .seealso: TSSetExactFinalTime() 1546 @*/ 1547 PetscErrorCode TSSetDuration(TS ts,PetscInt maxsteps,PetscReal maxtime) 1548 { 1549 PetscFunctionBegin; 1550 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1551 PetscValidLogicalCollectiveInt(ts,maxsteps,2); 1552 PetscValidLogicalCollectiveReal(ts,maxtime,2); 1553 if (maxsteps >= 0) ts->max_steps = maxsteps; 1554 if (maxtime != PETSC_DEFAULT) ts->max_time = maxtime; 1555 PetscFunctionReturn(0); 1556 } 1557 1558 #undef __FUNCT__ 1559 #define __FUNCT__ "TSSetSolution" 1560 /*@ 1561 TSSetSolution - Sets the initial solution vector 1562 for use by the TS routines. 1563 1564 Logically Collective on TS and Vec 1565 1566 Input Parameters: 1567 + ts - the TS context obtained from TSCreate() 1568 - x - the solution vector 1569 1570 Level: beginner 1571 1572 .keywords: TS, timestep, set, solution, initial conditions 1573 @*/ 1574 PetscErrorCode TSSetSolution(TS ts,Vec x) 1575 { 1576 PetscErrorCode ierr; 1577 DM dm; 1578 1579 PetscFunctionBegin; 1580 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1581 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 1582 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 1583 ierr = VecDestroy(&ts->vec_sol);CHKERRQ(ierr); 1584 ts->vec_sol = x; 1585 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1586 ierr = DMShellSetGlobalVector(dm,x);CHKERRQ(ierr); 1587 PetscFunctionReturn(0); 1588 } 1589 1590 #undef __FUNCT__ 1591 #define __FUNCT__ "TSSetPreStep" 1592 /*@C 1593 TSSetPreStep - Sets the general-purpose function 1594 called once at the beginning of each time step. 1595 1596 Logically Collective on TS 1597 1598 Input Parameters: 1599 + ts - The TS context obtained from TSCreate() 1600 - func - The function 1601 1602 Calling sequence of func: 1603 . func (TS ts); 1604 1605 Level: intermediate 1606 1607 Note: 1608 If a step is rejected, TSStep() will call this routine again before each attempt. 1609 The last completed time step number can be queried using TSGetTimeStepNumber(), the 1610 size of the step being attempted can be obtained using TSGetTimeStep(). 1611 1612 .keywords: TS, timestep 1613 .seealso: TSSetPreStage(), TSSetPostStep(), TSStep() 1614 @*/ 1615 PetscErrorCode TSSetPreStep(TS ts, PetscErrorCode (*func)(TS)) 1616 { 1617 PetscFunctionBegin; 1618 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1619 ts->ops->prestep = func; 1620 PetscFunctionReturn(0); 1621 } 1622 1623 #undef __FUNCT__ 1624 #define __FUNCT__ "TSPreStep" 1625 /*@ 1626 TSPreStep - Runs the user-defined pre-step function. 1627 1628 Collective on TS 1629 1630 Input Parameters: 1631 . ts - The TS context obtained from TSCreate() 1632 1633 Notes: 1634 TSPreStep() is typically used within time stepping implementations, 1635 so most users would not generally call this routine themselves. 1636 1637 Level: developer 1638 1639 .keywords: TS, timestep 1640 .seealso: TSSetPreStep(), TSPreStage(), TSPostStep() 1641 @*/ 1642 PetscErrorCode TSPreStep(TS ts) 1643 { 1644 PetscErrorCode ierr; 1645 1646 PetscFunctionBegin; 1647 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1648 if (ts->ops->prestep) { 1649 PetscStackPush("TS PreStep function"); 1650 ierr = (*ts->ops->prestep)(ts);CHKERRQ(ierr); 1651 PetscStackPop; 1652 } 1653 PetscFunctionReturn(0); 1654 } 1655 1656 #undef __FUNCT__ 1657 #define __FUNCT__ "TSSetPreStage" 1658 /*@C 1659 TSSetPreStage - Sets the general-purpose function 1660 called once at the beginning of each stage. 1661 1662 Logically Collective on TS 1663 1664 Input Parameters: 1665 + ts - The TS context obtained from TSCreate() 1666 - func - The function 1667 1668 Calling sequence of func: 1669 . PetscErrorCode func(TS ts, PetscReal stagetime); 1670 1671 Level: intermediate 1672 1673 Note: 1674 There may be several stages per time step. If the solve for a given stage fails, the step may be rejected and retried. 1675 The time step number being computed can be queried using TSGetTimeStepNumber() and the total size of the step being 1676 attempted can be obtained using TSGetTimeStep(). The time at the start of the step is available via TSGetTime(). 1677 1678 .keywords: TS, timestep 1679 .seealso: TSSetPreStep(), TSSetPostStep(), TSGetApplicationContext() 1680 @*/ 1681 PetscErrorCode TSSetPreStage(TS ts, PetscErrorCode (*func)(TS,PetscReal)) 1682 { 1683 PetscFunctionBegin; 1684 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1685 ts->ops->prestage = func; 1686 PetscFunctionReturn(0); 1687 } 1688 1689 #undef __FUNCT__ 1690 #define __FUNCT__ "TSPreStage" 1691 /*@ 1692 TSPreStage - Runs the user-defined pre-stage function set using TSSetPreStage() 1693 1694 Collective on TS 1695 1696 Input Parameters: 1697 . ts - The TS context obtained from TSCreate() 1698 1699 Notes: 1700 TSPreStage() is typically used within time stepping implementations, 1701 most users would not generally call this routine themselves. 1702 1703 Level: developer 1704 1705 .keywords: TS, timestep 1706 .seealso: TSSetPreStep(), TSPreStep(), TSPostStep() 1707 @*/ 1708 PetscErrorCode TSPreStage(TS ts, PetscReal stagetime) 1709 { 1710 PetscErrorCode ierr; 1711 1712 PetscFunctionBegin; 1713 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1714 if (ts->ops->prestage) { 1715 PetscStackPush("TS PreStage function"); 1716 ierr = (*ts->ops->prestage)(ts,stagetime);CHKERRQ(ierr); 1717 PetscStackPop; 1718 } 1719 PetscFunctionReturn(0); 1720 } 1721 1722 #undef __FUNCT__ 1723 #define __FUNCT__ "TSSetPostStep" 1724 /*@C 1725 TSSetPostStep - Sets the general-purpose function 1726 called once at the end of each time step. 1727 1728 Logically Collective on TS 1729 1730 Input Parameters: 1731 + ts - The TS context obtained from TSCreate() 1732 - func - The function 1733 1734 Calling sequence of func: 1735 $ func (TS ts); 1736 1737 Level: intermediate 1738 1739 .keywords: TS, timestep 1740 .seealso: TSSetPreStep(), TSSetPreStage(), TSGetTimeStep(), TSGetTimeStepNumber(), TSGetTime() 1741 @*/ 1742 PetscErrorCode TSSetPostStep(TS ts, PetscErrorCode (*func)(TS)) 1743 { 1744 PetscFunctionBegin; 1745 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1746 ts->ops->poststep = func; 1747 PetscFunctionReturn(0); 1748 } 1749 1750 #undef __FUNCT__ 1751 #define __FUNCT__ "TSPostStep" 1752 /*@ 1753 TSPostStep - Runs the user-defined post-step function. 1754 1755 Collective on TS 1756 1757 Input Parameters: 1758 . ts - The TS context obtained from TSCreate() 1759 1760 Notes: 1761 TSPostStep() is typically used within time stepping implementations, 1762 so most users would not generally call this routine themselves. 1763 1764 Level: developer 1765 1766 .keywords: TS, timestep 1767 @*/ 1768 PetscErrorCode TSPostStep(TS ts) 1769 { 1770 PetscErrorCode ierr; 1771 1772 PetscFunctionBegin; 1773 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1774 if (ts->ops->poststep) { 1775 PetscStackPush("TS PostStep function"); 1776 ierr = (*ts->ops->poststep)(ts);CHKERRQ(ierr); 1777 PetscStackPop; 1778 } 1779 PetscFunctionReturn(0); 1780 } 1781 1782 /* ------------ Routines to set performance monitoring options ----------- */ 1783 1784 #undef __FUNCT__ 1785 #define __FUNCT__ "TSMonitorSet" 1786 /*@C 1787 TSMonitorSet - Sets an ADDITIONAL function that is to be used at every 1788 timestep to display the iteration's progress. 1789 1790 Logically Collective on TS 1791 1792 Input Parameters: 1793 + ts - the TS context obtained from TSCreate() 1794 . monitor - monitoring routine 1795 . mctx - [optional] user-defined context for private data for the 1796 monitor routine (use PETSC_NULL if no context is desired) 1797 - monitordestroy - [optional] routine that frees monitor context 1798 (may be PETSC_NULL) 1799 1800 Calling sequence of monitor: 1801 $ int monitor(TS ts,PetscInt steps,PetscReal time,Vec x,void *mctx) 1802 1803 + ts - the TS context 1804 . steps - iteration number 1805 . time - current time 1806 . x - current iterate 1807 - mctx - [optional] monitoring context 1808 1809 Notes: 1810 This routine adds an additional monitor to the list of monitors that 1811 already has been loaded. 1812 1813 Fortran notes: Only a single monitor function can be set for each TS object 1814 1815 Level: intermediate 1816 1817 .keywords: TS, timestep, set, monitor 1818 1819 .seealso: TSMonitorDefault(), TSMonitorCancel() 1820 @*/ 1821 PetscErrorCode TSMonitorSet(TS ts,PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,void*),void *mctx,PetscErrorCode (*mdestroy)(void**)) 1822 { 1823 PetscFunctionBegin; 1824 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1825 if (ts->numbermonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1826 ts->monitor[ts->numbermonitors] = monitor; 1827 ts->mdestroy[ts->numbermonitors] = mdestroy; 1828 ts->monitorcontext[ts->numbermonitors++] = (void*)mctx; 1829 PetscFunctionReturn(0); 1830 } 1831 1832 #undef __FUNCT__ 1833 #define __FUNCT__ "TSMonitorCancel" 1834 /*@C 1835 TSMonitorCancel - Clears all the monitors that have been set on a time-step object. 1836 1837 Logically Collective on TS 1838 1839 Input Parameters: 1840 . ts - the TS context obtained from TSCreate() 1841 1842 Notes: 1843 There is no way to remove a single, specific monitor. 1844 1845 Level: intermediate 1846 1847 .keywords: TS, timestep, set, monitor 1848 1849 .seealso: TSMonitorDefault(), TSMonitorSet() 1850 @*/ 1851 PetscErrorCode TSMonitorCancel(TS ts) 1852 { 1853 PetscErrorCode ierr; 1854 PetscInt i; 1855 1856 PetscFunctionBegin; 1857 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1858 for (i=0; i<ts->numbermonitors; i++) { 1859 if (ts->mdestroy[i]) { 1860 ierr = (*ts->mdestroy[i])(&ts->monitorcontext[i]);CHKERRQ(ierr); 1861 } 1862 } 1863 ts->numbermonitors = 0; 1864 PetscFunctionReturn(0); 1865 } 1866 1867 #undef __FUNCT__ 1868 #define __FUNCT__ "TSMonitorDefault" 1869 /*@ 1870 TSMonitorDefault - Sets the Default monitor 1871 1872 Level: intermediate 1873 1874 .keywords: TS, set, monitor 1875 1876 .seealso: TSMonitorDefault(), TSMonitorSet() 1877 @*/ 1878 PetscErrorCode TSMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,void *dummy) 1879 { 1880 PetscErrorCode ierr; 1881 PetscViewer viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ts)->comm); 1882 1883 PetscFunctionBegin; 1884 ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1885 ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g\n",step,(double)ts->time_step,(double)ptime);CHKERRQ(ierr); 1886 ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1887 PetscFunctionReturn(0); 1888 } 1889 1890 #undef __FUNCT__ 1891 #define __FUNCT__ "TSSetRetainStages" 1892 /*@ 1893 TSSetRetainStages - Request that all stages in the upcoming step be stored so that interpolation will be available. 1894 1895 Logically Collective on TS 1896 1897 Input Argument: 1898 . ts - time stepping context 1899 1900 Output Argument: 1901 . flg - PETSC_TRUE or PETSC_FALSE 1902 1903 Level: intermediate 1904 1905 .keywords: TS, set 1906 1907 .seealso: TSInterpolate(), TSSetPostStep() 1908 @*/ 1909 PetscErrorCode TSSetRetainStages(TS ts,PetscBool flg) 1910 { 1911 1912 PetscFunctionBegin; 1913 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1914 ts->retain_stages = flg; 1915 PetscFunctionReturn(0); 1916 } 1917 1918 #undef __FUNCT__ 1919 #define __FUNCT__ "TSInterpolate" 1920 /*@ 1921 TSInterpolate - Interpolate the solution computed during the previous step to an arbitrary location in the interval 1922 1923 Collective on TS 1924 1925 Input Argument: 1926 + ts - time stepping context 1927 - t - time to interpolate to 1928 1929 Output Argument: 1930 . X - state at given time 1931 1932 Notes: 1933 The user should call TSSetRetainStages() before taking a step in which interpolation will be requested. 1934 1935 Level: intermediate 1936 1937 Developer Notes: 1938 TSInterpolate() and the storing of previous steps/stages should be generalized to support delay differential equations and continuous adjoints. 1939 1940 .keywords: TS, set 1941 1942 .seealso: TSSetRetainStages(), TSSetPostStep() 1943 @*/ 1944 PetscErrorCode TSInterpolate(TS ts,PetscReal t,Vec X) 1945 { 1946 PetscErrorCode ierr; 1947 1948 PetscFunctionBegin; 1949 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1950 if (t < ts->ptime - ts->time_step_prev || t > ts->ptime) SETERRQ3(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested time %G not in last time steps [%G,%G]",t,ts->ptime-ts->time_step_prev,ts->ptime); 1951 if (!ts->ops->interpolate) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"%s does not provide interpolation",((PetscObject)ts)->type_name); 1952 ierr = (*ts->ops->interpolate)(ts,t,X);CHKERRQ(ierr); 1953 PetscFunctionReturn(0); 1954 } 1955 1956 #undef __FUNCT__ 1957 #define __FUNCT__ "TSStep" 1958 /*@ 1959 TSStep - Steps one time step 1960 1961 Collective on TS 1962 1963 Input Parameter: 1964 . ts - the TS context obtained from TSCreate() 1965 1966 Level: intermediate 1967 1968 Notes: 1969 The hook set using TSSetPreStep() is called before each attempt to take the step. In general, the time step size may 1970 be changed due to adaptive error controller or solve failures. Note that steps may contain multiple stages. 1971 1972 .keywords: TS, timestep, solve 1973 1974 .seealso: TSCreate(), TSSetUp(), TSDestroy(), TSSolve(), TSSetPreStep(), TSSetPreStage() 1975 @*/ 1976 PetscErrorCode TSStep(TS ts) 1977 { 1978 PetscReal ptime_prev; 1979 PetscErrorCode ierr; 1980 1981 PetscFunctionBegin; 1982 PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1983 ierr = TSSetUp(ts);CHKERRQ(ierr); 1984 1985 ts->reason = TS_CONVERGED_ITERATING; 1986 1987 ptime_prev = ts->ptime; 1988 ierr = PetscLogEventBegin(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1989 ierr = (*ts->ops->step)(ts);CHKERRQ(ierr); 1990 ierr = PetscLogEventEnd(TS_Step,ts,0,0,0);CHKERRQ(ierr); 1991 ts->time_step_prev = ts->ptime - ptime_prev; 1992 1993 if (ts->reason < 0) { 1994 if (ts->errorifstepfailed) { 1995 if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) { 1996 SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]); 1997 } else SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]); 1998 } 1999 } else if (!ts->reason) { 2000 if (ts->steps >= ts->max_steps) 2001 ts->reason = TS_CONVERGED_ITS; 2002 else if (ts->ptime >= ts->max_time) 2003 ts->reason = TS_CONVERGED_TIME; 2004 } 2005 2006 PetscFunctionReturn(0); 2007 } 2008 2009 #undef __FUNCT__ 2010 #define __FUNCT__ "TSEvaluateStep" 2011 /*@ 2012 TSEvaluateStep - Evaluate the solution at the end of a time step with a given order of accuracy. 2013 2014 Collective on TS 2015 2016 Input Arguments: 2017 + ts - time stepping context 2018 . order - desired order of accuracy 2019 - done - whether the step was evaluated at this order (pass PETSC_NULL to generate an error if not available) 2020 2021 Output Arguments: 2022 . X - state at the end of the current step 2023 2024 Level: advanced 2025 2026 Notes: 2027 This function cannot be called until all stages have been evaluated. 2028 It is normally called by adaptive controllers before a step has been accepted and may also be called by the user after TSStep() has returned. 2029 2030 .seealso: TSStep(), TSAdapt 2031 @*/ 2032 PetscErrorCode TSEvaluateStep(TS ts,PetscInt order,Vec X,PetscBool *done) 2033 { 2034 PetscErrorCode ierr; 2035 2036 PetscFunctionBegin; 2037 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2038 PetscValidType(ts,1); 2039 PetscValidHeaderSpecific(X,VEC_CLASSID,3); 2040 if (!ts->ops->evaluatestep) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_SUP,"TSEvaluateStep not implemented for type '%s'",((PetscObject)ts)->type_name); 2041 ierr = (*ts->ops->evaluatestep)(ts,order,X,done);CHKERRQ(ierr); 2042 PetscFunctionReturn(0); 2043 } 2044 2045 #undef __FUNCT__ 2046 #define __FUNCT__ "TSSolve" 2047 /*@ 2048 TSSolve - Steps the requested number of timesteps. 2049 2050 Collective on TS 2051 2052 Input Parameter: 2053 + ts - the TS context obtained from TSCreate() 2054 - x - the solution vector 2055 2056 Output Parameter: 2057 . ftime - time of the state vector x upon completion 2058 2059 Level: beginner 2060 2061 Notes: 2062 The final time returned by this function may be different from the time of the internally 2063 held state accessible by TSGetSolution() and TSGetTime() because the method may have 2064 stepped over the final time. 2065 2066 .keywords: TS, timestep, solve 2067 2068 .seealso: TSCreate(), TSSetSolution(), TSStep() 2069 @*/ 2070 PetscErrorCode TSSolve(TS ts,Vec x,PetscReal *ftime) 2071 { 2072 PetscBool flg; 2073 char filename[PETSC_MAX_PATH_LEN]; 2074 PetscViewer viewer; 2075 PetscErrorCode ierr; 2076 2077 PetscFunctionBegin; 2078 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2079 PetscValidHeaderSpecific(x,VEC_CLASSID,2); 2080 if (ts->exact_final_time) { /* Need ts->vec_sol to be distinct so it is not overwritten when we interpolate at the end */ 2081 if (!ts->vec_sol || x == ts->vec_sol) { 2082 Vec y; 2083 ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 2084 ierr = TSSetSolution(ts,y);CHKERRQ(ierr); 2085 ierr = VecDestroy(&y);CHKERRQ(ierr); /* grant ownership */ 2086 } 2087 ierr = VecCopy(x,ts->vec_sol);CHKERRQ(ierr); 2088 } else { 2089 ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 2090 } 2091 ierr = TSSetUp(ts);CHKERRQ(ierr); 2092 /* reset time step and iteration counters */ 2093 ts->steps = 0; 2094 ts->ksp_its = 0; 2095 ts->snes_its = 0; 2096 ts->num_snes_failures = 0; 2097 ts->reject = 0; 2098 ts->reason = TS_CONVERGED_ITERATING; 2099 2100 if (ts->ops->solve) { /* This private interface is transitional and should be removed when all implementations are updated. */ 2101 ierr = (*ts->ops->solve)(ts);CHKERRQ(ierr); 2102 ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr); 2103 if (ftime) *ftime = ts->ptime; 2104 } else { 2105 /* steps the requested number of timesteps. */ 2106 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 2107 if (ts->steps >= ts->max_steps) 2108 ts->reason = TS_CONVERGED_ITS; 2109 else if (ts->ptime >= ts->max_time) 2110 ts->reason = TS_CONVERGED_TIME; 2111 while (!ts->reason) { 2112 ierr = TSStep(ts);CHKERRQ(ierr); 2113 ierr = TSPostStep(ts);CHKERRQ(ierr); 2114 ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 2115 } 2116 if (ts->exact_final_time && ts->ptime > ts->max_time) { 2117 ierr = TSInterpolate(ts,ts->max_time,x);CHKERRQ(ierr); 2118 if (ftime) *ftime = ts->max_time; 2119 } else { 2120 ierr = VecCopy(ts->vec_sol,x);CHKERRQ(ierr); 2121 if (ftime) *ftime = ts->ptime; 2122 } 2123 } 2124 ierr = PetscOptionsGetString(((PetscObject)ts)->prefix,"-ts_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2125 if (flg && !PetscPreLoadingOn) { 2126 ierr = PetscViewerASCIIOpen(((PetscObject)ts)->comm,filename,&viewer);CHKERRQ(ierr); 2127 ierr = TSView(ts,viewer);CHKERRQ(ierr); 2128 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 2129 } 2130 PetscFunctionReturn(0); 2131 } 2132 2133 #undef __FUNCT__ 2134 #define __FUNCT__ "TSMonitor" 2135 /*@ 2136 TSMonitor - Runs all user-provided monitor routines set using TSMonitorSet() 2137 2138 Collective on TS 2139 2140 Input Parameters: 2141 + ts - time stepping context obtained from TSCreate() 2142 . step - step number that has just completed 2143 . ptime - model time of the state 2144 - x - state at the current model time 2145 2146 Notes: 2147 TSMonitor() is typically used within the time stepping implementations. 2148 Users might call this function when using the TSStep() interface instead of TSSolve(). 2149 2150 Level: advanced 2151 2152 .keywords: TS, timestep 2153 @*/ 2154 PetscErrorCode TSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec x) 2155 { 2156 PetscErrorCode ierr; 2157 PetscInt i,n = ts->numbermonitors; 2158 2159 PetscFunctionBegin; 2160 for (i=0; i<n; i++) { 2161 ierr = (*ts->monitor[i])(ts,step,ptime,x,ts->monitorcontext[i]);CHKERRQ(ierr); 2162 } 2163 PetscFunctionReturn(0); 2164 } 2165 2166 /* ------------------------------------------------------------------------*/ 2167 2168 #undef __FUNCT__ 2169 #define __FUNCT__ "TSMonitorLGCreate" 2170 /*@C 2171 TSMonitorLGCreate - Creates a line graph context for use with 2172 TS to monitor convergence of preconditioned residual norms. 2173 2174 Collective on TS 2175 2176 Input Parameters: 2177 + host - the X display to open, or null for the local machine 2178 . label - the title to put in the title bar 2179 . x, y - the screen coordinates of the upper left coordinate of the window 2180 - m, n - the screen width and height in pixels 2181 2182 Output Parameter: 2183 . draw - the drawing context 2184 2185 Options Database Key: 2186 . -ts_monitor_draw - automatically sets line graph monitor 2187 2188 Notes: 2189 Use TSMonitorLGDestroy() to destroy this line graph, not PetscDrawLGDestroy(). 2190 2191 Level: intermediate 2192 2193 .keywords: TS, monitor, line graph, residual, seealso 2194 2195 .seealso: TSMonitorLGDestroy(), TSMonitorSet() 2196 2197 @*/ 2198 PetscErrorCode TSMonitorLGCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 2199 { 2200 PetscDraw win; 2201 PetscErrorCode ierr; 2202 PetscDrawAxis axis; 2203 2204 PetscFunctionBegin; 2205 ierr = PetscDrawCreate(comm,host,label,x,y,m,n,&win);CHKERRQ(ierr); 2206 ierr = PetscDrawSetType(win,PETSC_DRAW_X);CHKERRQ(ierr); 2207 ierr = PetscDrawLGCreate(win,1,draw);CHKERRQ(ierr); 2208 ierr = PetscDrawLGIndicateDataPoints(*draw);CHKERRQ(ierr); 2209 ierr = PetscDrawLGGetAxis(*draw,&axis);CHKERRQ(ierr); 2210 ierr = PetscDrawAxisSetLabels(axis,"Integration Time","Iteration","Time");CHKERRQ(ierr); 2211 ierr = PetscLogObjectParent(*draw,win);CHKERRQ(ierr); 2212 PetscFunctionReturn(0); 2213 } 2214 2215 #undef __FUNCT__ 2216 #define __FUNCT__ "TSMonitorLG" 2217 PetscErrorCode TSMonitorLG(TS ts,PetscInt n,PetscReal ptime,Vec v,void *monctx) 2218 { 2219 PetscDrawLG lg = (PetscDrawLG) monctx; 2220 PetscReal x,y = ptime; 2221 PetscErrorCode ierr; 2222 2223 PetscFunctionBegin; 2224 if (!monctx) { 2225 MPI_Comm comm; 2226 PetscViewer viewer; 2227 PetscDrawAxis axis; 2228 2229 ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr); 2230 viewer = PETSC_VIEWER_DRAW_(comm); 2231 ierr = PetscViewerDrawGetDrawLG(viewer,0,&lg);CHKERRQ(ierr); 2232 ierr = PetscDrawLGGetAxis(lg,&axis);CHKERRQ(ierr); 2233 ierr = PetscDrawAxisSetLabels(axis,"Integration Time","Iteration","Time");CHKERRQ(ierr); 2234 } 2235 2236 if (!n) {ierr = PetscDrawLGReset(lg);CHKERRQ(ierr);} 2237 x = (PetscReal)n; 2238 ierr = PetscDrawLGAddPoint(lg,&x,&y);CHKERRQ(ierr); 2239 if (n < 20 || (n % 5)) { 2240 ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 2241 } 2242 PetscFunctionReturn(0); 2243 } 2244 2245 #undef __FUNCT__ 2246 #define __FUNCT__ "TSMonitorLGDestroy" 2247 /*@C 2248 TSMonitorLGDestroy - Destroys a line graph context that was created 2249 with TSMonitorLGCreate(). 2250 2251 Collective on PetscDrawLG 2252 2253 Input Parameter: 2254 . draw - the drawing context 2255 2256 Level: intermediate 2257 2258 .keywords: TS, monitor, line graph, destroy 2259 2260 .seealso: TSMonitorLGCreate(), TSMonitorSet(), TSMonitorLG(); 2261 @*/ 2262 PetscErrorCode TSMonitorLGDestroy(PetscDrawLG *drawlg) 2263 { 2264 PetscDraw draw; 2265 PetscErrorCode ierr; 2266 2267 PetscFunctionBegin; 2268 ierr = PetscDrawLGGetDraw(*drawlg,&draw);CHKERRQ(ierr); 2269 ierr = PetscDrawDestroy(&draw);CHKERRQ(ierr); 2270 ierr = PetscDrawLGDestroy(drawlg);CHKERRQ(ierr); 2271 PetscFunctionReturn(0); 2272 } 2273 2274 #undef __FUNCT__ 2275 #define __FUNCT__ "TSGetTime" 2276 /*@ 2277 TSGetTime - Gets the time of the most recently completed step. 2278 2279 Not Collective 2280 2281 Input Parameter: 2282 . ts - the TS context obtained from TSCreate() 2283 2284 Output Parameter: 2285 . t - the current time 2286 2287 Level: beginner 2288 2289 Note: 2290 When called during time step evaluation (e.g. during residual evaluation or via hooks set using TSSetPreStep(), 2291 TSSetPreStage(), or TSSetPostStep()), the time is the time at the start of the step being evaluated. 2292 2293 .seealso: TSSetInitialTimeStep(), TSGetTimeStep() 2294 2295 .keywords: TS, get, time 2296 @*/ 2297 PetscErrorCode TSGetTime(TS ts,PetscReal* t) 2298 { 2299 PetscFunctionBegin; 2300 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2301 PetscValidRealPointer(t,2); 2302 *t = ts->ptime; 2303 PetscFunctionReturn(0); 2304 } 2305 2306 #undef __FUNCT__ 2307 #define __FUNCT__ "TSSetTime" 2308 /*@ 2309 TSSetTime - Allows one to reset the time. 2310 2311 Logically Collective on TS 2312 2313 Input Parameters: 2314 + ts - the TS context obtained from TSCreate() 2315 - time - the time 2316 2317 Level: intermediate 2318 2319 .seealso: TSGetTime(), TSSetDuration() 2320 2321 .keywords: TS, set, time 2322 @*/ 2323 PetscErrorCode TSSetTime(TS ts, PetscReal t) 2324 { 2325 PetscFunctionBegin; 2326 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2327 PetscValidLogicalCollectiveReal(ts,t,2); 2328 ts->ptime = t; 2329 PetscFunctionReturn(0); 2330 } 2331 2332 #undef __FUNCT__ 2333 #define __FUNCT__ "TSSetOptionsPrefix" 2334 /*@C 2335 TSSetOptionsPrefix - Sets the prefix used for searching for all 2336 TS options in the database. 2337 2338 Logically Collective on TS 2339 2340 Input Parameter: 2341 + ts - The TS context 2342 - prefix - The prefix to prepend to all option names 2343 2344 Notes: 2345 A hyphen (-) must NOT be given at the beginning of the prefix name. 2346 The first character of all runtime options is AUTOMATICALLY the 2347 hyphen. 2348 2349 Level: advanced 2350 2351 .keywords: TS, set, options, prefix, database 2352 2353 .seealso: TSSetFromOptions() 2354 2355 @*/ 2356 PetscErrorCode TSSetOptionsPrefix(TS ts,const char prefix[]) 2357 { 2358 PetscErrorCode ierr; 2359 SNES snes; 2360 2361 PetscFunctionBegin; 2362 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2363 ierr = PetscObjectSetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2364 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2365 ierr = SNESSetOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2366 PetscFunctionReturn(0); 2367 } 2368 2369 2370 #undef __FUNCT__ 2371 #define __FUNCT__ "TSAppendOptionsPrefix" 2372 /*@C 2373 TSAppendOptionsPrefix - Appends to the prefix used for searching for all 2374 TS options in the database. 2375 2376 Logically Collective on TS 2377 2378 Input Parameter: 2379 + ts - The TS context 2380 - prefix - The prefix to prepend to all option names 2381 2382 Notes: 2383 A hyphen (-) must NOT be given at the beginning of the prefix name. 2384 The first character of all runtime options is AUTOMATICALLY the 2385 hyphen. 2386 2387 Level: advanced 2388 2389 .keywords: TS, append, options, prefix, database 2390 2391 .seealso: TSGetOptionsPrefix() 2392 2393 @*/ 2394 PetscErrorCode TSAppendOptionsPrefix(TS ts,const char prefix[]) 2395 { 2396 PetscErrorCode ierr; 2397 SNES snes; 2398 2399 PetscFunctionBegin; 2400 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2401 ierr = PetscObjectAppendOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2402 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2403 ierr = SNESAppendOptionsPrefix(snes,prefix);CHKERRQ(ierr); 2404 PetscFunctionReturn(0); 2405 } 2406 2407 #undef __FUNCT__ 2408 #define __FUNCT__ "TSGetOptionsPrefix" 2409 /*@C 2410 TSGetOptionsPrefix - Sets the prefix used for searching for all 2411 TS options in the database. 2412 2413 Not Collective 2414 2415 Input Parameter: 2416 . ts - The TS context 2417 2418 Output Parameter: 2419 . prefix - A pointer to the prefix string used 2420 2421 Notes: On the fortran side, the user should pass in a string 'prifix' of 2422 sufficient length to hold the prefix. 2423 2424 Level: intermediate 2425 2426 .keywords: TS, get, options, prefix, database 2427 2428 .seealso: TSAppendOptionsPrefix() 2429 @*/ 2430 PetscErrorCode TSGetOptionsPrefix(TS ts,const char *prefix[]) 2431 { 2432 PetscErrorCode ierr; 2433 2434 PetscFunctionBegin; 2435 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2436 PetscValidPointer(prefix,2); 2437 ierr = PetscObjectGetOptionsPrefix((PetscObject)ts,prefix);CHKERRQ(ierr); 2438 PetscFunctionReturn(0); 2439 } 2440 2441 #undef __FUNCT__ 2442 #define __FUNCT__ "TSGetRHSJacobian" 2443 /*@C 2444 TSGetRHSJacobian - Returns the Jacobian J at the present timestep. 2445 2446 Not Collective, but parallel objects are returned if TS is parallel 2447 2448 Input Parameter: 2449 . ts - The TS context obtained from TSCreate() 2450 2451 Output Parameters: 2452 + J - The Jacobian J of F, where U_t = F(U,t) 2453 . M - The preconditioner matrix, usually the same as J 2454 . func - Function to compute the Jacobian of the RHS 2455 - ctx - User-defined context for Jacobian evaluation routine 2456 2457 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2458 2459 Level: intermediate 2460 2461 .seealso: TSGetTimeStep(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2462 2463 .keywords: TS, timestep, get, matrix, Jacobian 2464 @*/ 2465 PetscErrorCode TSGetRHSJacobian(TS ts,Mat *J,Mat *M,TSRHSJacobian *func,void **ctx) 2466 { 2467 PetscErrorCode ierr; 2468 SNES snes; 2469 DM dm; 2470 2471 PetscFunctionBegin; 2472 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2473 ierr = SNESGetJacobian(snes,J,M,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2474 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2475 ierr = DMTSGetRHSJacobian(dm,func,ctx);CHKERRQ(ierr); 2476 PetscFunctionReturn(0); 2477 } 2478 2479 #undef __FUNCT__ 2480 #define __FUNCT__ "TSGetIJacobian" 2481 /*@C 2482 TSGetIJacobian - Returns the implicit Jacobian at the present timestep. 2483 2484 Not Collective, but parallel objects are returned if TS is parallel 2485 2486 Input Parameter: 2487 . ts - The TS context obtained from TSCreate() 2488 2489 Output Parameters: 2490 + A - The Jacobian of F(t,U,U_t) 2491 . B - The preconditioner matrix, often the same as A 2492 . f - The function to compute the matrices 2493 - ctx - User-defined context for Jacobian evaluation routine 2494 2495 Notes: You can pass in PETSC_NULL for any return argument you do not need. 2496 2497 Level: advanced 2498 2499 .seealso: TSGetTimeStep(), TSGetRHSJacobian(), TSGetMatrices(), TSGetTime(), TSGetTimeStepNumber() 2500 2501 .keywords: TS, timestep, get, matrix, Jacobian 2502 @*/ 2503 PetscErrorCode TSGetIJacobian(TS ts,Mat *A,Mat *B,TSIJacobian *f,void **ctx) 2504 { 2505 PetscErrorCode ierr; 2506 SNES snes; 2507 DM dm; 2508 PetscFunctionBegin; 2509 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2510 ierr = SNESGetJacobian(snes,A,B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2511 ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 2512 ierr = DMTSGetIJacobian(dm,f,ctx);CHKERRQ(ierr); 2513 PetscFunctionReturn(0); 2514 } 2515 2516 typedef struct { 2517 PetscViewer viewer; 2518 Vec initialsolution; 2519 PetscBool showinitial; 2520 } TSMonitorSolutionCtx; 2521 2522 #undef __FUNCT__ 2523 #define __FUNCT__ "TSMonitorSolution" 2524 /*@C 2525 TSMonitorSolution - Monitors progress of the TS solvers by calling 2526 VecView() for the solution at each timestep 2527 2528 Collective on TS 2529 2530 Input Parameters: 2531 + ts - the TS context 2532 . step - current time-step 2533 . ptime - current time 2534 - dummy - either a viewer or PETSC_NULL 2535 2536 Level: intermediate 2537 2538 .keywords: TS, vector, monitor, view 2539 2540 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 2541 @*/ 2542 PetscErrorCode TSMonitorSolution(TS ts,PetscInt step,PetscReal ptime,Vec x,void *dummy) 2543 { 2544 PetscErrorCode ierr; 2545 TSMonitorSolutionCtx *ictx = (TSMonitorSolutionCtx*)dummy; 2546 2547 PetscFunctionBegin; 2548 if (!step && ictx->showinitial) { 2549 if (!ictx->initialsolution) { 2550 ierr = VecDuplicate(x,&ictx->initialsolution);CHKERRQ(ierr); 2551 } 2552 ierr = VecCopy(x,ictx->initialsolution);CHKERRQ(ierr); 2553 } 2554 if (ictx->showinitial) { 2555 PetscReal pause; 2556 ierr = PetscViewerDrawGetPause(ictx->viewer,&pause);CHKERRQ(ierr); 2557 ierr = PetscViewerDrawSetPause(ictx->viewer,0.0);CHKERRQ(ierr); 2558 ierr = VecView(ictx->initialsolution,ictx->viewer);CHKERRQ(ierr); 2559 ierr = PetscViewerDrawSetPause(ictx->viewer,pause);CHKERRQ(ierr); 2560 ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_TRUE);CHKERRQ(ierr); 2561 } 2562 ierr = VecView(x,ictx->viewer);CHKERRQ(ierr); 2563 if (ictx->showinitial) { 2564 ierr = PetscViewerDrawSetHold(ictx->viewer,PETSC_FALSE);CHKERRQ(ierr); 2565 } 2566 PetscFunctionReturn(0); 2567 } 2568 2569 2570 #undef __FUNCT__ 2571 #define __FUNCT__ "TSMonitorSolutionDestroy" 2572 /*@C 2573 TSMonitorSolutionDestroy - Destroys the monitor context for TSMonitorSolution 2574 2575 Collective on TS 2576 2577 Input Parameters: 2578 . ctx - the monitor context 2579 2580 Level: intermediate 2581 2582 .keywords: TS, vector, monitor, view 2583 2584 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution() 2585 @*/ 2586 PetscErrorCode TSMonitorSolutionDestroy(void **ctx) 2587 { 2588 PetscErrorCode ierr; 2589 TSMonitorSolutionCtx *ictx = *(TSMonitorSolutionCtx**)ctx; 2590 2591 PetscFunctionBegin; 2592 ierr = PetscViewerDestroy(&ictx->viewer);CHKERRQ(ierr); 2593 ierr = VecDestroy(&ictx->initialsolution);CHKERRQ(ierr); 2594 ierr = PetscFree(ictx);CHKERRQ(ierr); 2595 PetscFunctionReturn(0); 2596 } 2597 2598 #undef __FUNCT__ 2599 #define __FUNCT__ "TSMonitorSolutionCreate" 2600 /*@C 2601 TSMonitorSolutionCreate - Creates the monitor context for TSMonitorSolution 2602 2603 Collective on TS 2604 2605 Input Parameter: 2606 . ts - time-step context 2607 2608 Output Patameter: 2609 . ctx - the monitor context 2610 2611 Level: intermediate 2612 2613 .keywords: TS, vector, monitor, view 2614 2615 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView(), TSMonitorSolution() 2616 @*/ 2617 PetscErrorCode TSMonitorSolutionCreate(TS ts,PetscViewer viewer,void **ctx) 2618 { 2619 PetscErrorCode ierr; 2620 TSMonitorSolutionCtx *ictx; 2621 2622 PetscFunctionBegin; 2623 ierr = PetscNew(TSMonitorSolutionCtx,&ictx);CHKERRQ(ierr); 2624 *ctx = (void*)ictx; 2625 if (!viewer) { 2626 viewer = PETSC_VIEWER_DRAW_(((PetscObject)ts)->comm); 2627 } 2628 ierr = PetscObjectReference((PetscObject)viewer);CHKERRQ(ierr); 2629 ictx->viewer = viewer; 2630 ictx->showinitial = PETSC_FALSE; 2631 ierr = PetscOptionsGetBool(((PetscObject)ts)->prefix,"-ts_monitor_solution_initial",&ictx->showinitial,PETSC_NULL);CHKERRQ(ierr); 2632 PetscFunctionReturn(0); 2633 } 2634 2635 #undef __FUNCT__ 2636 #define __FUNCT__ "TSSetDM" 2637 /*@ 2638 TSSetDM - Sets the DM that may be used by some preconditioners 2639 2640 Logically Collective on TS and DM 2641 2642 Input Parameters: 2643 + ts - the preconditioner context 2644 - dm - the dm 2645 2646 Level: intermediate 2647 2648 2649 .seealso: TSGetDM(), SNESSetDM(), SNESGetDM() 2650 @*/ 2651 PetscErrorCode TSSetDM(TS ts,DM dm) 2652 { 2653 PetscErrorCode ierr; 2654 SNES snes; 2655 TSDM tsdm; 2656 2657 PetscFunctionBegin; 2658 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2659 ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr); 2660 if (ts->dm) { /* Move the TSDM context over to the new DM unless the new DM already has one */ 2661 PetscContainer oldcontainer,container; 2662 ierr = PetscObjectQuery((PetscObject)ts->dm,"TSDM",(PetscObject*)&oldcontainer);CHKERRQ(ierr); 2663 ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr); 2664 if (oldcontainer && !container) { 2665 ierr = DMTSCopyContext(ts->dm,dm);CHKERRQ(ierr); 2666 ierr = DMTSGetContext(ts->dm,&tsdm);CHKERRQ(ierr); 2667 if (tsdm->originaldm == ts->dm) { /* Grant write privileges to the replacement DM */ 2668 tsdm->originaldm = dm; 2669 } 2670 } 2671 ierr = DMDestroy(&ts->dm);CHKERRQ(ierr); 2672 } 2673 ts->dm = dm; 2674 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 2675 ierr = SNESSetDM(snes,dm);CHKERRQ(ierr); 2676 PetscFunctionReturn(0); 2677 } 2678 2679 #undef __FUNCT__ 2680 #define __FUNCT__ "TSGetDM" 2681 /*@ 2682 TSGetDM - Gets the DM that may be used by some preconditioners 2683 2684 Not Collective 2685 2686 Input Parameter: 2687 . ts - the preconditioner context 2688 2689 Output Parameter: 2690 . dm - the dm 2691 2692 Level: intermediate 2693 2694 2695 .seealso: TSSetDM(), SNESSetDM(), SNESGetDM() 2696 @*/ 2697 PetscErrorCode TSGetDM(TS ts,DM *dm) 2698 { 2699 PetscErrorCode ierr; 2700 2701 PetscFunctionBegin; 2702 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2703 if (!ts->dm) { 2704 ierr = DMShellCreate(((PetscObject)ts)->comm,&ts->dm);CHKERRQ(ierr); 2705 if (ts->snes) {ierr = SNESSetDM(ts->snes,ts->dm);CHKERRQ(ierr);} 2706 } 2707 *dm = ts->dm; 2708 PetscFunctionReturn(0); 2709 } 2710 2711 #undef __FUNCT__ 2712 #define __FUNCT__ "SNESTSFormFunction" 2713 /*@ 2714 SNESTSFormFunction - Function to evaluate nonlinear residual 2715 2716 Logically Collective on SNES 2717 2718 Input Parameter: 2719 + snes - nonlinear solver 2720 . X - the current state at which to evaluate the residual 2721 - ctx - user context, must be a TS 2722 2723 Output Parameter: 2724 . F - the nonlinear residual 2725 2726 Notes: 2727 This function is not normally called by users and is automatically registered with the SNES used by TS. 2728 It is most frequently passed to MatFDColoringSetFunction(). 2729 2730 Level: advanced 2731 2732 .seealso: SNESSetFunction(), MatFDColoringSetFunction() 2733 @*/ 2734 PetscErrorCode SNESTSFormFunction(SNES snes,Vec X,Vec F,void *ctx) 2735 { 2736 TS ts = (TS)ctx; 2737 PetscErrorCode ierr; 2738 2739 PetscFunctionBegin; 2740 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2741 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2742 PetscValidHeaderSpecific(F,VEC_CLASSID,3); 2743 PetscValidHeaderSpecific(ts,TS_CLASSID,4); 2744 ierr = (ts->ops->snesfunction)(snes,X,F,ts);CHKERRQ(ierr); 2745 PetscFunctionReturn(0); 2746 } 2747 2748 #undef __FUNCT__ 2749 #define __FUNCT__ "SNESTSFormJacobian" 2750 /*@ 2751 SNESTSFormJacobian - Function to evaluate the Jacobian 2752 2753 Collective on SNES 2754 2755 Input Parameter: 2756 + snes - nonlinear solver 2757 . X - the current state at which to evaluate the residual 2758 - ctx - user context, must be a TS 2759 2760 Output Parameter: 2761 + A - the Jacobian 2762 . B - the preconditioning matrix (may be the same as A) 2763 - flag - indicates any structure change in the matrix 2764 2765 Notes: 2766 This function is not normally called by users and is automatically registered with the SNES used by TS. 2767 2768 Level: developer 2769 2770 .seealso: SNESSetJacobian() 2771 @*/ 2772 PetscErrorCode SNESTSFormJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flag,void *ctx) 2773 { 2774 TS ts = (TS)ctx; 2775 PetscErrorCode ierr; 2776 2777 PetscFunctionBegin; 2778 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 2779 PetscValidHeaderSpecific(X,VEC_CLASSID,2); 2780 PetscValidPointer(A,3); 2781 PetscValidHeaderSpecific(*A,MAT_CLASSID,3); 2782 PetscValidPointer(B,4); 2783 PetscValidHeaderSpecific(*B,MAT_CLASSID,4); 2784 PetscValidPointer(flag,5); 2785 PetscValidHeaderSpecific(ts,TS_CLASSID,6); 2786 ierr = (ts->ops->snesjacobian)(snes,X,A,B,flag,ts);CHKERRQ(ierr); 2787 PetscFunctionReturn(0); 2788 } 2789 2790 #undef __FUNCT__ 2791 #define __FUNCT__ "TSComputeRHSFunctionLinear" 2792 /*@C 2793 TSComputeRHSFunctionLinear - Evaluate the right hand side via the user-provided Jacobian, for linear problems only 2794 2795 Collective on TS 2796 2797 Input Arguments: 2798 + ts - time stepping context 2799 . t - time at which to evaluate 2800 . X - state at which to evaluate 2801 - ctx - context 2802 2803 Output Arguments: 2804 . F - right hand side 2805 2806 Level: intermediate 2807 2808 Notes: 2809 This function is intended to be passed to TSSetRHSFunction() to evaluate the right hand side for linear problems. 2810 The matrix (and optionally the evaluation context) should be passed to TSSetRHSJacobian(). 2811 2812 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSJacobianConstant() 2813 @*/ 2814 PetscErrorCode TSComputeRHSFunctionLinear(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 2815 { 2816 PetscErrorCode ierr; 2817 Mat Arhs,Brhs; 2818 MatStructure flg2; 2819 2820 PetscFunctionBegin; 2821 ierr = TSGetRHSMats_Private(ts,&Arhs,&Brhs);CHKERRQ(ierr); 2822 ierr = TSComputeRHSJacobian(ts,t,X,&Arhs,&Brhs,&flg2);CHKERRQ(ierr); 2823 ierr = MatMult(Arhs,X,F);CHKERRQ(ierr); 2824 PetscFunctionReturn(0); 2825 } 2826 2827 #undef __FUNCT__ 2828 #define __FUNCT__ "TSComputeRHSJacobianConstant" 2829 /*@C 2830 TSComputeRHSJacobianConstant - Reuses a Jacobian that is time-independent. 2831 2832 Collective on TS 2833 2834 Input Arguments: 2835 + ts - time stepping context 2836 . t - time at which to evaluate 2837 . X - state at which to evaluate 2838 - ctx - context 2839 2840 Output Arguments: 2841 + A - pointer to operator 2842 . B - pointer to preconditioning matrix 2843 - flg - matrix structure flag 2844 2845 Level: intermediate 2846 2847 Notes: 2848 This function is intended to be passed to TSSetRHSJacobian() to evaluate the Jacobian for linear time-independent problems. 2849 2850 .seealso: TSSetRHSFunction(), TSSetRHSJacobian(), TSComputeRHSFunctionLinear() 2851 @*/ 2852 PetscErrorCode TSComputeRHSJacobianConstant(TS ts,PetscReal t,Vec X,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2853 { 2854 2855 PetscFunctionBegin; 2856 *flg = SAME_PRECONDITIONER; 2857 PetscFunctionReturn(0); 2858 } 2859 2860 #undef __FUNCT__ 2861 #define __FUNCT__ "TSComputeIFunctionLinear" 2862 /*@C 2863 TSComputeIFunctionLinear - Evaluate the left hand side via the user-provided Jacobian, for linear problems only 2864 2865 Collective on TS 2866 2867 Input Arguments: 2868 + ts - time stepping context 2869 . t - time at which to evaluate 2870 . X - state at which to evaluate 2871 . Xdot - time derivative of state vector 2872 - ctx - context 2873 2874 Output Arguments: 2875 . F - left hand side 2876 2877 Level: intermediate 2878 2879 Notes: 2880 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 2881 user is required to write their own TSComputeIFunction. 2882 This function is intended to be passed to TSSetIFunction() to evaluate the left hand side for linear problems. 2883 The matrix (and optionally the evaluation context) should be passed to TSSetIJacobian(). 2884 2885 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIJacobianConstant() 2886 @*/ 2887 PetscErrorCode TSComputeIFunctionLinear(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 2888 { 2889 PetscErrorCode ierr; 2890 Mat A,B; 2891 MatStructure flg2; 2892 2893 PetscFunctionBegin; 2894 ierr = TSGetIJacobian(ts,&A,&B,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 2895 ierr = TSComputeIJacobian(ts,t,X,Xdot,1.0,&A,&B,&flg2,PETSC_TRUE);CHKERRQ(ierr); 2896 ierr = MatMult(A,Xdot,F);CHKERRQ(ierr); 2897 PetscFunctionReturn(0); 2898 } 2899 2900 #undef __FUNCT__ 2901 #define __FUNCT__ "TSComputeIJacobianConstant" 2902 /*@C 2903 TSComputeIJacobianConstant - Reuses a Jacobian that is time-independent. 2904 2905 Collective on TS 2906 2907 Input Arguments: 2908 + ts - time stepping context 2909 . t - time at which to evaluate 2910 . X - state at which to evaluate 2911 . Xdot - time derivative of state vector 2912 . shift - shift to apply 2913 - ctx - context 2914 2915 Output Arguments: 2916 + A - pointer to operator 2917 . B - pointer to preconditioning matrix 2918 - flg - matrix structure flag 2919 2920 Level: intermediate 2921 2922 Notes: 2923 This function is intended to be passed to TSSetIJacobian() to evaluate the Jacobian for linear time-independent problems. 2924 2925 .seealso: TSSetIFunction(), TSSetIJacobian(), TSComputeIFunctionLinear() 2926 @*/ 2927 PetscErrorCode TSComputeIJacobianConstant(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flg,void *ctx) 2928 { 2929 2930 PetscFunctionBegin; 2931 *flg = SAME_PRECONDITIONER; 2932 PetscFunctionReturn(0); 2933 } 2934 2935 2936 #undef __FUNCT__ 2937 #define __FUNCT__ "TSGetConvergedReason" 2938 /*@ 2939 TSGetConvergedReason - Gets the reason the TS iteration was stopped. 2940 2941 Not Collective 2942 2943 Input Parameter: 2944 . ts - the TS context 2945 2946 Output Parameter: 2947 . reason - negative value indicates diverged, positive value converged, see TSConvergedReason or the 2948 manual pages for the individual convergence tests for complete lists 2949 2950 Level: intermediate 2951 2952 Notes: 2953 Can only be called after the call to TSSolve() is complete. 2954 2955 .keywords: TS, nonlinear, set, convergence, test 2956 2957 .seealso: TSSetConvergenceTest(), TSConvergedReason 2958 @*/ 2959 PetscErrorCode TSGetConvergedReason(TS ts,TSConvergedReason *reason) 2960 { 2961 PetscFunctionBegin; 2962 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2963 PetscValidPointer(reason,2); 2964 *reason = ts->reason; 2965 PetscFunctionReturn(0); 2966 } 2967 2968 #undef __FUNCT__ 2969 #define __FUNCT__ "TSGetSNESIterations" 2970 /*@ 2971 TSGetSNESIterations - Gets the total number of nonlinear iterations 2972 used by the time integrator. 2973 2974 Not Collective 2975 2976 Input Parameter: 2977 . ts - TS context 2978 2979 Output Parameter: 2980 . nits - number of nonlinear iterations 2981 2982 Notes: 2983 This counter is reset to zero for each successive call to TSSolve(). 2984 2985 Level: intermediate 2986 2987 .keywords: TS, get, number, nonlinear, iterations 2988 2989 .seealso: TSGetKSPIterations() 2990 @*/ 2991 PetscErrorCode TSGetSNESIterations(TS ts,PetscInt *nits) 2992 { 2993 PetscFunctionBegin; 2994 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2995 PetscValidIntPointer(nits,2); 2996 *nits = ts->snes_its; 2997 PetscFunctionReturn(0); 2998 } 2999 3000 #undef __FUNCT__ 3001 #define __FUNCT__ "TSGetKSPIterations" 3002 /*@ 3003 TSGetKSPIterations - Gets the total number of linear iterations 3004 used by the time integrator. 3005 3006 Not Collective 3007 3008 Input Parameter: 3009 . ts - TS context 3010 3011 Output Parameter: 3012 . lits - number of linear iterations 3013 3014 Notes: 3015 This counter is reset to zero for each successive call to TSSolve(). 3016 3017 Level: intermediate 3018 3019 .keywords: TS, get, number, linear, iterations 3020 3021 .seealso: TSGetSNESIterations(), SNESGetKSPIterations() 3022 @*/ 3023 PetscErrorCode TSGetKSPIterations(TS ts,PetscInt *lits) 3024 { 3025 PetscFunctionBegin; 3026 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3027 PetscValidIntPointer(lits,2); 3028 *lits = ts->ksp_its; 3029 PetscFunctionReturn(0); 3030 } 3031 3032 #undef __FUNCT__ 3033 #define __FUNCT__ "TSGetStepRejections" 3034 /*@ 3035 TSGetStepRejections - Gets the total number of rejected steps. 3036 3037 Not Collective 3038 3039 Input Parameter: 3040 . ts - TS context 3041 3042 Output Parameter: 3043 . rejects - number of steps rejected 3044 3045 Notes: 3046 This counter is reset to zero for each successive call to TSSolve(). 3047 3048 Level: intermediate 3049 3050 .keywords: TS, get, number 3051 3052 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetSNESFailures(), TSSetMaxSNESFailures(), TSSetErrorIfStepFails() 3053 @*/ 3054 PetscErrorCode TSGetStepRejections(TS ts,PetscInt *rejects) 3055 { 3056 PetscFunctionBegin; 3057 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3058 PetscValidIntPointer(rejects,2); 3059 *rejects = ts->reject; 3060 PetscFunctionReturn(0); 3061 } 3062 3063 #undef __FUNCT__ 3064 #define __FUNCT__ "TSGetSNESFailures" 3065 /*@ 3066 TSGetSNESFailures - Gets the total number of failed SNES solves 3067 3068 Not Collective 3069 3070 Input Parameter: 3071 . ts - TS context 3072 3073 Output Parameter: 3074 . fails - number of failed nonlinear solves 3075 3076 Notes: 3077 This counter is reset to zero for each successive call to TSSolve(). 3078 3079 Level: intermediate 3080 3081 .keywords: TS, get, number 3082 3083 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSSetMaxSNESFailures() 3084 @*/ 3085 PetscErrorCode TSGetSNESFailures(TS ts,PetscInt *fails) 3086 { 3087 PetscFunctionBegin; 3088 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3089 PetscValidIntPointer(fails,2); 3090 *fails = ts->num_snes_failures; 3091 PetscFunctionReturn(0); 3092 } 3093 3094 #undef __FUNCT__ 3095 #define __FUNCT__ "TSSetMaxStepRejections" 3096 /*@ 3097 TSSetMaxStepRejections - Sets the maximum number of step rejections before a step fails 3098 3099 Not Collective 3100 3101 Input Parameter: 3102 + ts - TS context 3103 - rejects - maximum number of rejected steps, pass -1 for unlimited 3104 3105 Notes: 3106 The counter is reset to zero for each step 3107 3108 Options Database Key: 3109 . -ts_max_reject - Maximum number of step rejections before a step fails 3110 3111 Level: intermediate 3112 3113 .keywords: TS, set, maximum, number 3114 3115 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxSNESFailures(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason() 3116 @*/ 3117 PetscErrorCode TSSetMaxStepRejections(TS ts,PetscInt rejects) 3118 { 3119 PetscFunctionBegin; 3120 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3121 ts->max_reject = rejects; 3122 PetscFunctionReturn(0); 3123 } 3124 3125 #undef __FUNCT__ 3126 #define __FUNCT__ "TSSetMaxSNESFailures" 3127 /*@ 3128 TSSetMaxSNESFailures - Sets the maximum number of failed SNES solves 3129 3130 Not Collective 3131 3132 Input Parameter: 3133 + ts - TS context 3134 - fails - maximum number of failed nonlinear solves, pass -1 for unlimited 3135 3136 Notes: 3137 The counter is reset to zero for each successive call to TSSolve(). 3138 3139 Options Database Key: 3140 . -ts_max_snes_failures - Maximum number of nonlinear solve failures 3141 3142 Level: intermediate 3143 3144 .keywords: TS, set, maximum, number 3145 3146 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), SNESGetConvergedReason(), TSGetConvergedReason() 3147 @*/ 3148 PetscErrorCode TSSetMaxSNESFailures(TS ts,PetscInt fails) 3149 { 3150 PetscFunctionBegin; 3151 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3152 ts->max_snes_failures = fails; 3153 PetscFunctionReturn(0); 3154 } 3155 3156 #undef __FUNCT__ 3157 #define __FUNCT__ "TSSetErrorIfStepFails()" 3158 /*@ 3159 TSSetErrorIfStepFails - Error if no step succeeds 3160 3161 Not Collective 3162 3163 Input Parameter: 3164 + ts - TS context 3165 - err - PETSC_TRUE to error if no step succeeds, PETSC_FALSE to return without failure 3166 3167 Options Database Key: 3168 . -ts_error_if_step_fails - Error if no step succeeds 3169 3170 Level: intermediate 3171 3172 .keywords: TS, set, error 3173 3174 .seealso: TSGetSNESIterations(), TSGetKSPIterations(), TSSetMaxStepRejections(), TSGetStepRejections(), TSGetSNESFailures(), TSSetErrorIfStepFails(), TSGetConvergedReason() 3175 @*/ 3176 PetscErrorCode TSSetErrorIfStepFails(TS ts,PetscBool err) 3177 { 3178 PetscFunctionBegin; 3179 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3180 ts->errorifstepfailed = err; 3181 PetscFunctionReturn(0); 3182 } 3183 3184 #undef __FUNCT__ 3185 #define __FUNCT__ "TSMonitorSolutionBinary" 3186 /*@C 3187 TSMonitorSolutionBinary - Monitors progress of the TS solvers by VecView() for the solution at each timestep. Normally the viewer is a binary file 3188 3189 Collective on TS 3190 3191 Input Parameters: 3192 + ts - the TS context 3193 . step - current time-step 3194 . ptime - current time 3195 . x - current state 3196 - viewer - binary viewer 3197 3198 Level: intermediate 3199 3200 .keywords: TS, vector, monitor, view 3201 3202 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 3203 @*/ 3204 PetscErrorCode TSMonitorSolutionBinary(TS ts,PetscInt step,PetscReal ptime,Vec x,void *viewer) 3205 { 3206 PetscErrorCode ierr; 3207 PetscViewer v = (PetscViewer)viewer; 3208 3209 PetscFunctionBegin; 3210 ierr = VecView(x,v);CHKERRQ(ierr); 3211 PetscFunctionReturn(0); 3212 } 3213 3214 #undef __FUNCT__ 3215 #define __FUNCT__ "TSMonitorSolutionVTK" 3216 /*@C 3217 TSMonitorSolutionVTK - Monitors progress of the TS solvers by VecView() for the solution at each timestep. 3218 3219 Collective on TS 3220 3221 Input Parameters: 3222 + ts - the TS context 3223 . step - current time-step 3224 . ptime - current time 3225 . x - current state 3226 - filenametemplate - string containing a format specifier for the integer time step (e.g. %03D) 3227 3228 Level: intermediate 3229 3230 Notes: 3231 The VTK format does not allow writing multiple time steps in the same file, therefore a different file will be written for each time step. 3232 These are named according to the file name template. 3233 3234 This function is normally passed as an argument to TSMonitorSet() along with TSMonitorSolutionVTKDestroy(). 3235 3236 .keywords: TS, vector, monitor, view 3237 3238 .seealso: TSMonitorSet(), TSMonitorDefault(), VecView() 3239 @*/ 3240 PetscErrorCode TSMonitorSolutionVTK(TS ts,PetscInt step,PetscReal ptime,Vec x,void *filenametemplate) 3241 { 3242 PetscErrorCode ierr; 3243 char filename[PETSC_MAX_PATH_LEN]; 3244 PetscViewer viewer; 3245 3246 PetscFunctionBegin; 3247 ierr = PetscSNPrintf(filename,sizeof(filename),(const char*)filenametemplate,step);CHKERRQ(ierr); 3248 ierr = PetscViewerVTKOpen(((PetscObject)ts)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 3249 ierr = VecView(x,viewer);CHKERRQ(ierr); 3250 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 3251 PetscFunctionReturn(0); 3252 } 3253 3254 #undef __FUNCT__ 3255 #define __FUNCT__ "TSMonitorSolutionVTKDestroy" 3256 /*@C 3257 TSMonitorSolutionVTKDestroy - Destroy context for monitoring 3258 3259 Collective on TS 3260 3261 Input Parameters: 3262 . filenametemplate - string containing a format specifier for the integer time step (e.g. %03D) 3263 3264 Level: intermediate 3265 3266 Note: 3267 This function is normally passed to TSMonitorSet() along with TSMonitorSolutionVTK(). 3268 3269 .keywords: TS, vector, monitor, view 3270 3271 .seealso: TSMonitorSet(), TSMonitorSolutionVTK() 3272 @*/ 3273 PetscErrorCode TSMonitorSolutionVTKDestroy(void *filenametemplate) 3274 { 3275 PetscErrorCode ierr; 3276 3277 PetscFunctionBegin; 3278 ierr = PetscFree(*(char**)filenametemplate);CHKERRQ(ierr); 3279 PetscFunctionReturn(0); 3280 } 3281 3282 #undef __FUNCT__ 3283 #define __FUNCT__ "TSGetAdapt" 3284 /*@ 3285 TSGetAdapt - Get the adaptive controller context for the current method 3286 3287 Collective on TS if controller has not been created yet 3288 3289 Input Arguments: 3290 . ts - time stepping context 3291 3292 Output Arguments: 3293 . adapt - adaptive controller 3294 3295 Level: intermediate 3296 3297 .seealso: TSAdapt, TSAdaptSetType(), TSAdaptChoose() 3298 @*/ 3299 PetscErrorCode TSGetAdapt(TS ts,TSAdapt *adapt) 3300 { 3301 PetscErrorCode ierr; 3302 3303 PetscFunctionBegin; 3304 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3305 PetscValidPointer(adapt,2); 3306 if (!ts->adapt) { 3307 ierr = TSAdaptCreate(((PetscObject)ts)->comm,&ts->adapt);CHKERRQ(ierr); 3308 ierr = PetscLogObjectParent(ts,ts->adapt);CHKERRQ(ierr); 3309 ierr = PetscObjectIncrementTabLevel((PetscObject)ts->adapt,(PetscObject)ts,1);CHKERRQ(ierr); 3310 } 3311 *adapt = ts->adapt; 3312 PetscFunctionReturn(0); 3313 } 3314 3315 #undef __FUNCT__ 3316 #define __FUNCT__ "TSSetTolerances" 3317 /*@ 3318 TSSetTolerances - Set tolerances for local truncation error when using adaptive controller 3319 3320 Logically Collective 3321 3322 Input Arguments: 3323 + ts - time integration context 3324 . atol - scalar absolute tolerances, PETSC_DECIDE to leave current value 3325 . vatol - vector of absolute tolerances or PETSC_NULL, used in preference to atol if present 3326 . rtol - scalar relative tolerances, PETSC_DECIDE to leave current value 3327 - vrtol - vector of relative tolerances or PETSC_NULL, used in preference to atol if present 3328 3329 Level: beginner 3330 3331 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSGetTolerances() 3332 @*/ 3333 PetscErrorCode TSSetTolerances(TS ts,PetscReal atol,Vec vatol,PetscReal rtol,Vec vrtol) 3334 { 3335 PetscErrorCode ierr; 3336 3337 PetscFunctionBegin; 3338 if (atol != PETSC_DECIDE && atol != PETSC_DEFAULT) ts->atol = atol; 3339 if (vatol) { 3340 ierr = PetscObjectReference((PetscObject)vatol);CHKERRQ(ierr); 3341 ierr = VecDestroy(&ts->vatol);CHKERRQ(ierr); 3342 ts->vatol = vatol; 3343 } 3344 if (rtol != PETSC_DECIDE && rtol != PETSC_DEFAULT) ts->rtol = rtol; 3345 if (vrtol) { 3346 ierr = PetscObjectReference((PetscObject)vrtol);CHKERRQ(ierr); 3347 ierr = VecDestroy(&ts->vrtol);CHKERRQ(ierr); 3348 ts->vrtol = vrtol; 3349 } 3350 PetscFunctionReturn(0); 3351 } 3352 3353 #undef __FUNCT__ 3354 #define __FUNCT__ "TSGetTolerances" 3355 /*@ 3356 TSGetTolerances - Get tolerances for local truncation error when using adaptive controller 3357 3358 Logically Collective 3359 3360 Input Arguments: 3361 . ts - time integration context 3362 3363 Output Arguments: 3364 + atol - scalar absolute tolerances, PETSC_NULL to ignore 3365 . vatol - vector of absolute tolerances, PETSC_NULL to ignore 3366 . rtol - scalar relative tolerances, PETSC_NULL to ignore 3367 - vrtol - vector of relative tolerances, PETSC_NULL to ignore 3368 3369 Level: beginner 3370 3371 .seealso: TS, TSAdapt, TSVecNormWRMS(), TSSetTolerances() 3372 @*/ 3373 PetscErrorCode TSGetTolerances(TS ts,PetscReal *atol,Vec *vatol,PetscReal *rtol,Vec *vrtol) 3374 { 3375 3376 PetscFunctionBegin; 3377 if (atol) *atol = ts->atol; 3378 if (vatol) *vatol = ts->vatol; 3379 if (rtol) *rtol = ts->rtol; 3380 if (vrtol) *vrtol = ts->vrtol; 3381 PetscFunctionReturn(0); 3382 } 3383 3384 #undef __FUNCT__ 3385 #define __FUNCT__ "TSErrorNormWRMS" 3386 /*@ 3387 TSErrorNormWRMS - compute a weighted norm of the difference between a vector and the current state 3388 3389 Collective on TS 3390 3391 Input Arguments: 3392 + ts - time stepping context 3393 - Y - state vector to be compared to ts->vec_sol 3394 3395 Output Arguments: 3396 . norm - weighted norm, a value of 1.0 is considered small 3397 3398 Level: developer 3399 3400 .seealso: TSSetTolerances() 3401 @*/ 3402 PetscErrorCode TSErrorNormWRMS(TS ts,Vec Y,PetscReal *norm) 3403 { 3404 PetscErrorCode ierr; 3405 PetscInt i,n,N; 3406 const PetscScalar *x,*y; 3407 Vec X; 3408 PetscReal sum,gsum; 3409 3410 PetscFunctionBegin; 3411 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3412 PetscValidHeaderSpecific(Y,VEC_CLASSID,2); 3413 PetscValidPointer(norm,3); 3414 X = ts->vec_sol; 3415 PetscCheckSameTypeAndComm(X,1,Y,2); 3416 if (X == Y) SETERRQ(((PetscObject)X)->comm,PETSC_ERR_ARG_IDN,"Y cannot be the TS solution vector"); 3417 3418 ierr = VecGetSize(X,&N);CHKERRQ(ierr); 3419 ierr = VecGetLocalSize(X,&n);CHKERRQ(ierr); 3420 ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 3421 ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 3422 sum = 0.; 3423 if (ts->vatol && ts->vrtol) { 3424 const PetscScalar *atol,*rtol; 3425 ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3426 ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3427 for (i=0; i<n; i++) { 3428 PetscReal tol = PetscRealPart(atol[i]) + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3429 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3430 } 3431 ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3432 ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3433 } else if (ts->vatol) { /* vector atol, scalar rtol */ 3434 const PetscScalar *atol; 3435 ierr = VecGetArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3436 for (i=0; i<n; i++) { 3437 PetscReal tol = PetscRealPart(atol[i]) + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3438 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3439 } 3440 ierr = VecRestoreArrayRead(ts->vatol,&atol);CHKERRQ(ierr); 3441 } else if (ts->vrtol) { /* scalar atol, vector rtol */ 3442 const PetscScalar *rtol; 3443 ierr = VecGetArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3444 for (i=0; i<n; i++) { 3445 PetscReal tol = ts->atol + PetscRealPart(rtol[i]) * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3446 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3447 } 3448 ierr = VecRestoreArrayRead(ts->vrtol,&rtol);CHKERRQ(ierr); 3449 } else { /* scalar atol, scalar rtol */ 3450 for (i=0; i<n; i++) { 3451 PetscReal tol = ts->atol + ts->rtol * PetscMax(PetscAbsScalar(x[i]),PetscAbsScalar(y[i])); 3452 sum += PetscSqr(PetscAbsScalar(y[i] - x[i]) / tol); 3453 } 3454 } 3455 ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 3456 ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 3457 3458 ierr = MPI_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,((PetscObject)ts)->comm);CHKERRQ(ierr); 3459 *norm = PetscSqrtReal(gsum / N); 3460 if (PetscIsInfOrNanScalar(*norm)) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 3461 PetscFunctionReturn(0); 3462 } 3463 3464 #undef __FUNCT__ 3465 #define __FUNCT__ "TSSetCFLTimeLocal" 3466 /*@ 3467 TSSetCFLTimeLocal - Set the local CFL constraint relative to forward Euler 3468 3469 Logically Collective on TS 3470 3471 Input Arguments: 3472 + ts - time stepping context 3473 - cfltime - maximum stable time step if using forward Euler (value can be different on each process) 3474 3475 Note: 3476 After calling this function, the global CFL time can be obtained by calling TSGetCFLTime() 3477 3478 Level: intermediate 3479 3480 .seealso: TSGetCFLTime(), TSADAPTCFL 3481 @*/ 3482 PetscErrorCode TSSetCFLTimeLocal(TS ts,PetscReal cfltime) 3483 { 3484 3485 PetscFunctionBegin; 3486 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3487 ts->cfltime_local = cfltime; 3488 ts->cfltime = -1.; 3489 PetscFunctionReturn(0); 3490 } 3491 3492 #undef __FUNCT__ 3493 #define __FUNCT__ "TSGetCFLTime" 3494 /*@ 3495 TSGetCFLTime - Get the maximum stable time step according to CFL criteria applied to forward Euler 3496 3497 Collective on TS 3498 3499 Input Arguments: 3500 . ts - time stepping context 3501 3502 Output Arguments: 3503 . cfltime - maximum stable time step for forward Euler 3504 3505 Level: advanced 3506 3507 .seealso: TSSetCFLTimeLocal() 3508 @*/ 3509 PetscErrorCode TSGetCFLTime(TS ts,PetscReal *cfltime) 3510 { 3511 PetscErrorCode ierr; 3512 3513 PetscFunctionBegin; 3514 if (ts->cfltime < 0) { 3515 ierr = MPI_Allreduce(&ts->cfltime_local,&ts->cfltime,1,MPIU_REAL,MPIU_MIN,((PetscObject)ts)->comm);CHKERRQ(ierr); 3516 } 3517 *cfltime = ts->cfltime; 3518 PetscFunctionReturn(0); 3519 } 3520 3521 #undef __FUNCT__ 3522 #define __FUNCT__ "TSVISetVariableBounds" 3523 /*@ 3524 TSVISetVariableBounds - Sets the lower and upper bounds for the solution vector. xl <= x <= xu 3525 3526 Input Parameters: 3527 . ts - the TS context. 3528 . xl - lower bound. 3529 . xu - upper bound. 3530 3531 Notes: 3532 If this routine is not called then the lower and upper bounds are set to 3533 SNES_VI_NINF and SNES_VI_INF respectively during SNESSetUp(). 3534 3535 Level: advanced 3536 3537 @*/ 3538 PetscErrorCode TSVISetVariableBounds(TS ts, Vec xl, Vec xu) 3539 { 3540 PetscErrorCode ierr; 3541 SNES snes; 3542 3543 PetscFunctionBegin; 3544 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 3545 ierr = SNESVISetVariableBounds(snes,xl,xu);CHKERRQ(ierr); 3546 PetscFunctionReturn(0); 3547 } 3548 3549 #if defined(PETSC_HAVE_MATLAB_ENGINE) 3550 #include <mex.h> 3551 3552 typedef struct {char *funcname; mxArray *ctx;} TSMatlabContext; 3553 3554 #undef __FUNCT__ 3555 #define __FUNCT__ "TSComputeFunction_Matlab" 3556 /* 3557 TSComputeFunction_Matlab - Calls the function that has been set with 3558 TSSetFunctionMatlab(). 3559 3560 Collective on TS 3561 3562 Input Parameters: 3563 + snes - the TS context 3564 - x - input vector 3565 3566 Output Parameter: 3567 . y - function vector, as set by TSSetFunction() 3568 3569 Notes: 3570 TSComputeFunction() is typically used within nonlinear solvers 3571 implementations, so most users would not generally call this routine 3572 themselves. 3573 3574 Level: developer 3575 3576 .keywords: TS, nonlinear, compute, function 3577 3578 .seealso: TSSetFunction(), TSGetFunction() 3579 */ 3580 PetscErrorCode TSComputeFunction_Matlab(TS snes,PetscReal time,Vec x,Vec xdot,Vec y, void *ctx) 3581 { 3582 PetscErrorCode ierr; 3583 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 3584 int nlhs = 1,nrhs = 7; 3585 mxArray *plhs[1],*prhs[7]; 3586 long long int lx = 0,lxdot = 0,ly = 0,ls = 0; 3587 3588 PetscFunctionBegin; 3589 PetscValidHeaderSpecific(snes,TS_CLASSID,1); 3590 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 3591 PetscValidHeaderSpecific(xdot,VEC_CLASSID,4); 3592 PetscValidHeaderSpecific(y,VEC_CLASSID,5); 3593 PetscCheckSameComm(snes,1,x,3); 3594 PetscCheckSameComm(snes,1,y,5); 3595 3596 ierr = PetscMemcpy(&ls,&snes,sizeof(snes));CHKERRQ(ierr); 3597 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3598 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(xdot));CHKERRQ(ierr); 3599 ierr = PetscMemcpy(&ly,&y,sizeof(x));CHKERRQ(ierr); 3600 prhs[0] = mxCreateDoubleScalar((double)ls); 3601 prhs[1] = mxCreateDoubleScalar(time); 3602 prhs[2] = mxCreateDoubleScalar((double)lx); 3603 prhs[3] = mxCreateDoubleScalar((double)lxdot); 3604 prhs[4] = mxCreateDoubleScalar((double)ly); 3605 prhs[5] = mxCreateString(sctx->funcname); 3606 prhs[6] = sctx->ctx; 3607 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeFunctionInternal");CHKERRQ(ierr); 3608 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3609 mxDestroyArray(prhs[0]); 3610 mxDestroyArray(prhs[1]); 3611 mxDestroyArray(prhs[2]); 3612 mxDestroyArray(prhs[3]); 3613 mxDestroyArray(prhs[4]); 3614 mxDestroyArray(prhs[5]); 3615 mxDestroyArray(plhs[0]); 3616 PetscFunctionReturn(0); 3617 } 3618 3619 3620 #undef __FUNCT__ 3621 #define __FUNCT__ "TSSetFunctionMatlab" 3622 /* 3623 TSSetFunctionMatlab - Sets the function evaluation routine and function 3624 vector for use by the TS routines in solving ODEs 3625 equations from MATLAB. Here the function is a string containing the name of a MATLAB function 3626 3627 Logically Collective on TS 3628 3629 Input Parameters: 3630 + ts - the TS context 3631 - func - function evaluation routine 3632 3633 Calling sequence of func: 3634 $ func (TS ts,PetscReal time,Vec x,Vec xdot,Vec f,void *ctx); 3635 3636 Level: beginner 3637 3638 .keywords: TS, nonlinear, set, function 3639 3640 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 3641 */ 3642 PetscErrorCode TSSetFunctionMatlab(TS ts,const char *func,mxArray *ctx) 3643 { 3644 PetscErrorCode ierr; 3645 TSMatlabContext *sctx; 3646 3647 PetscFunctionBegin; 3648 /* currently sctx is memory bleed */ 3649 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 3650 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3651 /* 3652 This should work, but it doesn't 3653 sctx->ctx = ctx; 3654 mexMakeArrayPersistent(sctx->ctx); 3655 */ 3656 sctx->ctx = mxDuplicateArray(ctx); 3657 ierr = TSSetIFunction(ts,PETSC_NULL,TSComputeFunction_Matlab,sctx);CHKERRQ(ierr); 3658 PetscFunctionReturn(0); 3659 } 3660 3661 #undef __FUNCT__ 3662 #define __FUNCT__ "TSComputeJacobian_Matlab" 3663 /* 3664 TSComputeJacobian_Matlab - Calls the function that has been set with 3665 TSSetJacobianMatlab(). 3666 3667 Collective on TS 3668 3669 Input Parameters: 3670 + ts - the TS context 3671 . x - input vector 3672 . A, B - the matrices 3673 - ctx - user context 3674 3675 Output Parameter: 3676 . flag - structure of the matrix 3677 3678 Level: developer 3679 3680 .keywords: TS, nonlinear, compute, function 3681 3682 .seealso: TSSetFunction(), TSGetFunction() 3683 @*/ 3684 PetscErrorCode TSComputeJacobian_Matlab(TS ts,PetscReal time,Vec x,Vec xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *flag, void *ctx) 3685 { 3686 PetscErrorCode ierr; 3687 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 3688 int nlhs = 2,nrhs = 9; 3689 mxArray *plhs[2],*prhs[9]; 3690 long long int lx = 0,lxdot = 0,lA = 0,ls = 0, lB = 0; 3691 3692 PetscFunctionBegin; 3693 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3694 PetscValidHeaderSpecific(x,VEC_CLASSID,3); 3695 3696 /* call Matlab function in ctx with arguments x and y */ 3697 3698 ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr); 3699 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3700 ierr = PetscMemcpy(&lxdot,&xdot,sizeof(x));CHKERRQ(ierr); 3701 ierr = PetscMemcpy(&lA,A,sizeof(x));CHKERRQ(ierr); 3702 ierr = PetscMemcpy(&lB,B,sizeof(x));CHKERRQ(ierr); 3703 prhs[0] = mxCreateDoubleScalar((double)ls); 3704 prhs[1] = mxCreateDoubleScalar((double)time); 3705 prhs[2] = mxCreateDoubleScalar((double)lx); 3706 prhs[3] = mxCreateDoubleScalar((double)lxdot); 3707 prhs[4] = mxCreateDoubleScalar((double)shift); 3708 prhs[5] = mxCreateDoubleScalar((double)lA); 3709 prhs[6] = mxCreateDoubleScalar((double)lB); 3710 prhs[7] = mxCreateString(sctx->funcname); 3711 prhs[8] = sctx->ctx; 3712 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSComputeJacobianInternal");CHKERRQ(ierr); 3713 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3714 *flag = (MatStructure) mxGetScalar(plhs[1]);CHKERRQ(ierr); 3715 mxDestroyArray(prhs[0]); 3716 mxDestroyArray(prhs[1]); 3717 mxDestroyArray(prhs[2]); 3718 mxDestroyArray(prhs[3]); 3719 mxDestroyArray(prhs[4]); 3720 mxDestroyArray(prhs[5]); 3721 mxDestroyArray(prhs[6]); 3722 mxDestroyArray(prhs[7]); 3723 mxDestroyArray(plhs[0]); 3724 mxDestroyArray(plhs[1]); 3725 PetscFunctionReturn(0); 3726 } 3727 3728 3729 #undef __FUNCT__ 3730 #define __FUNCT__ "TSSetJacobianMatlab" 3731 /* 3732 TSSetJacobianMatlab - Sets the Jacobian function evaluation routine and two empty Jacobian matrices 3733 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 3734 3735 Logically Collective on TS 3736 3737 Input Parameters: 3738 + ts - the TS context 3739 . A,B - Jacobian matrices 3740 . func - function evaluation routine 3741 - ctx - user context 3742 3743 Calling sequence of func: 3744 $ flag = func (TS ts,PetscReal time,Vec x,Vec xdot,Mat A,Mat B,void *ctx); 3745 3746 3747 Level: developer 3748 3749 .keywords: TS, nonlinear, set, function 3750 3751 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 3752 */ 3753 PetscErrorCode TSSetJacobianMatlab(TS ts,Mat A,Mat B,const char *func,mxArray *ctx) 3754 { 3755 PetscErrorCode ierr; 3756 TSMatlabContext *sctx; 3757 3758 PetscFunctionBegin; 3759 /* currently sctx is memory bleed */ 3760 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 3761 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3762 /* 3763 This should work, but it doesn't 3764 sctx->ctx = ctx; 3765 mexMakeArrayPersistent(sctx->ctx); 3766 */ 3767 sctx->ctx = mxDuplicateArray(ctx); 3768 ierr = TSSetIJacobian(ts,A,B,TSComputeJacobian_Matlab,sctx);CHKERRQ(ierr); 3769 PetscFunctionReturn(0); 3770 } 3771 3772 #undef __FUNCT__ 3773 #define __FUNCT__ "TSMonitor_Matlab" 3774 /* 3775 TSMonitor_Matlab - Calls the function that has been set with TSMonitorSetMatlab(). 3776 3777 Collective on TS 3778 3779 .seealso: TSSetFunction(), TSGetFunction() 3780 @*/ 3781 PetscErrorCode TSMonitor_Matlab(TS ts,PetscInt it, PetscReal time,Vec x, void *ctx) 3782 { 3783 PetscErrorCode ierr; 3784 TSMatlabContext *sctx = (TSMatlabContext *)ctx; 3785 int nlhs = 1,nrhs = 6; 3786 mxArray *plhs[1],*prhs[6]; 3787 long long int lx = 0,ls = 0; 3788 3789 PetscFunctionBegin; 3790 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 3791 PetscValidHeaderSpecific(x,VEC_CLASSID,4); 3792 3793 ierr = PetscMemcpy(&ls,&ts,sizeof(ts));CHKERRQ(ierr); 3794 ierr = PetscMemcpy(&lx,&x,sizeof(x));CHKERRQ(ierr); 3795 prhs[0] = mxCreateDoubleScalar((double)ls); 3796 prhs[1] = mxCreateDoubleScalar((double)it); 3797 prhs[2] = mxCreateDoubleScalar((double)time); 3798 prhs[3] = mxCreateDoubleScalar((double)lx); 3799 prhs[4] = mxCreateString(sctx->funcname); 3800 prhs[5] = sctx->ctx; 3801 ierr = mexCallMATLAB(nlhs,plhs,nrhs,prhs,"PetscTSMonitorInternal");CHKERRQ(ierr); 3802 ierr = mxGetScalar(plhs[0]);CHKERRQ(ierr); 3803 mxDestroyArray(prhs[0]); 3804 mxDestroyArray(prhs[1]); 3805 mxDestroyArray(prhs[2]); 3806 mxDestroyArray(prhs[3]); 3807 mxDestroyArray(prhs[4]); 3808 mxDestroyArray(plhs[0]); 3809 PetscFunctionReturn(0); 3810 } 3811 3812 3813 #undef __FUNCT__ 3814 #define __FUNCT__ "TSMonitorSetMatlab" 3815 /* 3816 TSMonitorSetMatlab - Sets the monitor function from Matlab 3817 3818 Level: developer 3819 3820 .keywords: TS, nonlinear, set, function 3821 3822 .seealso: TSGetFunction(), TSComputeFunction(), TSSetJacobian(), TSSetFunction() 3823 */ 3824 PetscErrorCode TSMonitorSetMatlab(TS ts,const char *func,mxArray *ctx) 3825 { 3826 PetscErrorCode ierr; 3827 TSMatlabContext *sctx; 3828 3829 PetscFunctionBegin; 3830 /* currently sctx is memory bleed */ 3831 ierr = PetscMalloc(sizeof(TSMatlabContext),&sctx);CHKERRQ(ierr); 3832 ierr = PetscStrallocpy(func,&sctx->funcname);CHKERRQ(ierr); 3833 /* 3834 This should work, but it doesn't 3835 sctx->ctx = ctx; 3836 mexMakeArrayPersistent(sctx->ctx); 3837 */ 3838 sctx->ctx = mxDuplicateArray(ctx); 3839 ierr = TSMonitorSet(ts,TSMonitor_Matlab,sctx,PETSC_NULL);CHKERRQ(ierr); 3840 PetscFunctionReturn(0); 3841 } 3842 #endif 3843