1 /* 2 User interface for the timestepping package. This package 3 is for use in solving time-dependent PDEs. 4 */ 5 #if !defined(__PETSCTS_H) 6 #define __PETSCTS_H 7 #include "petscsnes.h" 8 PETSC_EXTERN_CXX_BEGIN 9 10 /*S 11 TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 12 13 Level: beginner 14 15 Concepts: ODE solvers 16 17 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC 18 S*/ 19 typedef struct _p_TS* TS; 20 21 /*J 22 TSType - String with the name of a PETSc TS method or the creation function 23 with an optional dynamic library name, for example 24 http://www.mcs.anl.gov/petsc/lib.a:mytscreate() 25 26 Level: beginner 27 28 .seealso: TSSetType(), TS 29 J*/ 30 #define TSType char* 31 #define TSEULER "euler" 32 #define TSBEULER "beuler" 33 #define TSPSEUDO "pseudo" 34 #define TSCN "cn" 35 #define TSSUNDIALS "sundials" 36 #define TSRK "rk" 37 #define TSPYTHON "python" 38 #define TSTHETA "theta" 39 #define TSALPHA "alpha" 40 #define TSGL "gl" 41 #define TSSSP "ssp" 42 #define TSARKIMEX "arkimex" 43 #define TSROSW "rosw" 44 45 /*E 46 TSProblemType - Determines the type of problem this TS object is to be used to solve 47 48 Level: beginner 49 50 .seealso: TSCreate() 51 E*/ 52 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType; 53 54 /*E 55 TSConvergedReason - reason a TS method has converged or not 56 57 Level: beginner 58 59 Developer Notes: this must match finclude/petscts.h 60 61 Each reason has its own manual page. 62 63 .seealso: TSGetConvergedReason() 64 E*/ 65 typedef enum { 66 TS_CONVERGED_ITERATING = 0, 67 TS_CONVERGED_TIME = 1, 68 TS_CONVERGED_ITS = 2, 69 TS_DIVERGED_NONLINEAR_SOLVE = -1, 70 TS_DIVERGED_STEP_REJECTED = -2 71 } TSConvergedReason; 72 extern const char *const*TSConvergedReasons; 73 74 /*MC 75 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 76 77 Level: beginner 78 79 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt() 80 M*/ 81 82 /*MC 83 TS_CONVERGED_TIME - the final time was reached 84 85 Level: beginner 86 87 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration() 88 M*/ 89 90 /*MC 91 TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time 92 93 Level: beginner 94 95 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration() 96 M*/ 97 98 /*MC 99 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 100 101 Level: beginner 102 103 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason() 104 M*/ 105 106 /*MC 107 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 108 109 Level: beginner 110 111 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt() 112 M*/ 113 114 /* Logging support */ 115 extern PetscClassId TS_CLASSID; 116 117 extern PetscErrorCode TSInitializePackage(const char[]); 118 119 extern PetscErrorCode TSCreate(MPI_Comm,TS*); 120 extern PetscErrorCode TSDestroy(TS*); 121 122 extern PetscErrorCode TSSetProblemType(TS,TSProblemType); 123 extern PetscErrorCode TSGetProblemType(TS,TSProblemType*); 124 extern PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 125 extern PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 126 extern PetscErrorCode TSMonitorCancel(TS); 127 128 extern PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 129 extern PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 130 extern PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 131 extern PetscErrorCode TSSetFromOptions(TS); 132 extern PetscErrorCode TSSetUp(TS); 133 extern PetscErrorCode TSReset(TS); 134 135 extern PetscErrorCode TSSetSolution(TS,Vec); 136 extern PetscErrorCode TSGetSolution(TS,Vec*); 137 138 extern PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 139 extern PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 140 extern PetscErrorCode TSSetExactFinalTime(TS,PetscBool); 141 142 extern PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 143 extern PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*); 144 extern PetscErrorCode TSMonitorSolutionCreate(TS,PetscViewer,void**); 145 extern PetscErrorCode TSMonitorSolutionDestroy(void**); 146 extern PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 147 extern PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 148 extern PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 149 150 extern PetscErrorCode TSStep(TS); 151 extern PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 152 extern PetscErrorCode TSSolve(TS,Vec,PetscReal*); 153 extern PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 154 extern PetscErrorCode TSGetNonlinearSolveIterations(TS,PetscInt*); 155 extern PetscErrorCode TSGetLinearSolveIterations(TS,PetscInt*); 156 157 extern PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 158 extern PetscErrorCode TSGetTimeStep(TS,PetscReal*); 159 extern PetscErrorCode TSGetTime(TS,PetscReal*); 160 extern PetscErrorCode TSSetTime(TS,PetscReal); 161 extern PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 162 extern PetscErrorCode TSSetTimeStep(TS,PetscReal); 163 164 typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 165 typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 166 extern PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 167 extern PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 168 extern PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 169 extern PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 170 171 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 172 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 173 extern PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 174 extern PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 175 extern PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 176 extern PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 177 178 extern PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 179 extern PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 180 extern PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 181 extern PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 182 183 extern PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 184 extern PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 185 extern PetscErrorCode TSPreStep(TS); 186 extern PetscErrorCode TSPostStep(TS); 187 extern PetscErrorCode TSSetRetainStages(TS,PetscBool); 188 extern PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 189 extern PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 190 extern PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 191 extern PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 192 extern PetscErrorCode TSGetCFLTime(TS,PetscReal*); 193 194 extern PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 195 extern PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 196 extern PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 197 extern PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 198 199 extern PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 200 extern PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 201 extern PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 202 extern PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 203 extern PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 204 205 extern PetscErrorCode TSPythonSetType(TS,const char[]); 206 207 extern PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 208 extern PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 209 extern PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 210 extern PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 211 212 extern PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 213 214 /* Dynamic creation and loading functions */ 215 extern PetscFList TSList; 216 extern PetscBool TSRegisterAllCalled; 217 extern PetscErrorCode TSGetType(TS,const TSType*); 218 extern PetscErrorCode TSSetType(TS,const TSType); 219 extern PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 220 extern PetscErrorCode TSRegisterAll(const char[]); 221 extern PetscErrorCode TSRegisterDestroy(void); 222 223 /*MC 224 TSRegisterDynamic - Adds a creation method to the TS package. 225 226 Synopsis: 227 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 228 229 Not Collective 230 231 Input Parameters: 232 + name - The name of a new user-defined creation routine 233 . path - The path (either absolute or relative) of the library containing this routine 234 . func_name - The name of the creation routine 235 - create_func - The creation routine itself 236 237 Notes: 238 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 239 240 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 241 242 Sample usage: 243 .vb 244 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 245 .ve 246 247 Then, your ts type can be chosen with the procedural interface via 248 .vb 249 TS ts; 250 TSCreate(MPI_Comm, &ts); 251 TSSetType(ts, "my_ts") 252 .ve 253 or at runtime via the option 254 .vb 255 -ts_type my_ts 256 .ve 257 258 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 259 If your function is not being put into a shared library then use TSRegister() instead 260 261 Level: advanced 262 263 .keywords: TS, register 264 .seealso: TSRegisterAll(), TSRegisterDestroy() 265 M*/ 266 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 267 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 268 #else 269 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 270 #endif 271 272 extern PetscErrorCode TSGetSNES(TS,SNES*); 273 extern PetscErrorCode TSGetKSP(TS,KSP*); 274 275 extern PetscErrorCode TSView(TS,PetscViewer); 276 277 extern PetscErrorCode TSSetApplicationContext(TS,void *); 278 extern PetscErrorCode TSGetApplicationContext(TS,void *); 279 280 extern PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *); 281 extern PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *); 282 extern PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*); 283 284 /*J 285 TSSSPType - string with the name of TSSSP scheme. 286 287 Level: beginner 288 289 .seealso: TSSSPSetType(), TS 290 J*/ 291 #define TSSSPType char* 292 #define TSSSPRKS2 "rks2" 293 #define TSSSPRKS3 "rks3" 294 #define TSSSPRK104 "rk104" 295 296 extern PetscErrorCode TSSSPSetType(TS,const TSSSPType); 297 extern PetscErrorCode TSSSPGetType(TS,const TSSSPType*); 298 extern PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 299 extern PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 300 301 /*S 302 TSAdapt - Abstract object that manages time-step adaptivity 303 304 Level: beginner 305 306 .seealso: TS, TSAdaptCreate(), TSAdaptType 307 S*/ 308 typedef struct _p_TSAdapt *TSAdapt; 309 310 /*E 311 TSAdaptType - String with the name of TSAdapt scheme or the creation function 312 with an optional dynamic library name, for example 313 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 314 315 Level: beginner 316 317 .seealso: TSAdaptSetType(), TS 318 E*/ 319 #define TSAdaptType char* 320 #define TSADAPTBASIC "basic" 321 #define TSADAPTNONE "none" 322 #define TSADAPTCFL "cfl" 323 324 /*MC 325 TSAdaptRegisterDynamic - adds a TSAdapt implementation 326 327 Synopsis: 328 PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 329 330 Not Collective 331 332 Input Parameters: 333 + name_scheme - name of user-defined adaptivity scheme 334 . path - path (either absolute or relative) the library containing this scheme 335 . name_create - name of routine to create method context 336 - routine_create - routine to create method context 337 338 Notes: 339 TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 340 341 If dynamic libraries are used, then the fourth input argument (routine_create) 342 is ignored. 343 344 Sample usage: 345 .vb 346 TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 347 "MySchemeCreate",MySchemeCreate); 348 .ve 349 350 Then, your scheme can be chosen with the procedural interface via 351 $ TSAdaptSetType(ts,"my_scheme") 352 or at runtime via the option 353 $ -ts_adapt_type my_scheme 354 355 Level: advanced 356 357 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 358 and others of the form ${any_environmental_variable} occuring in pathname will be 359 replaced with appropriate values. 360 361 .keywords: TSAdapt, register 362 363 .seealso: TSAdaptRegisterAll() 364 M*/ 365 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 366 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0) 367 #else 368 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d) 369 #endif 370 371 extern PetscErrorCode TSGetAdapt(TS,TSAdapt*); 372 extern PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt)); 373 extern PetscErrorCode TSAdaptRegisterAll(const char[]); 374 extern PetscErrorCode TSAdaptRegisterDestroy(void); 375 extern PetscErrorCode TSAdaptInitializePackage(const char[]); 376 extern PetscErrorCode TSAdaptFinalizePackage(void); 377 extern PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 378 extern PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType); 379 extern PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 380 extern PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 381 extern PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 382 extern PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 383 extern PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 384 extern PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 385 extern PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 386 extern PetscErrorCode TSAdaptDestroy(TSAdapt*); 387 extern PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 388 extern PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 389 390 /*S 391 TSGLAdapt - Abstract object that manages time-step adaptivity 392 393 Level: beginner 394 395 Developer Notes: 396 This functionality should be replaced by the TSAdapt. 397 398 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 399 S*/ 400 typedef struct _p_TSGLAdapt *TSGLAdapt; 401 402 /*J 403 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 404 with an optional dynamic library name, for example 405 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 406 407 Level: beginner 408 409 .seealso: TSGLAdaptSetType(), TS 410 J*/ 411 #define TSGLAdaptType char* 412 #define TSGLADAPT_NONE "none" 413 #define TSGLADAPT_SIZE "size" 414 #define TSGLADAPT_BOTH "both" 415 416 /*MC 417 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 418 419 Synopsis: 420 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 421 422 Not Collective 423 424 Input Parameters: 425 + name_scheme - name of user-defined adaptivity scheme 426 . path - path (either absolute or relative) the library containing this scheme 427 . name_create - name of routine to create method context 428 - routine_create - routine to create method context 429 430 Notes: 431 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 432 433 If dynamic libraries are used, then the fourth input argument (routine_create) 434 is ignored. 435 436 Sample usage: 437 .vb 438 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 439 "MySchemeCreate",MySchemeCreate); 440 .ve 441 442 Then, your scheme can be chosen with the procedural interface via 443 $ TSGLAdaptSetType(ts,"my_scheme") 444 or at runtime via the option 445 $ -ts_adapt_type my_scheme 446 447 Level: advanced 448 449 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 450 and others of the form ${any_environmental_variable} occuring in pathname will be 451 replaced with appropriate values. 452 453 .keywords: TSGLAdapt, register 454 455 .seealso: TSGLAdaptRegisterAll() 456 M*/ 457 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 458 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 459 #else 460 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 461 #endif 462 463 extern PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 464 extern PetscErrorCode TSGLAdaptRegisterAll(const char[]); 465 extern PetscErrorCode TSGLAdaptRegisterDestroy(void); 466 extern PetscErrorCode TSGLAdaptInitializePackage(const char[]); 467 extern PetscErrorCode TSGLAdaptFinalizePackage(void); 468 extern PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 469 extern PetscErrorCode TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType); 470 extern PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 471 extern PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 472 extern PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 473 extern PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 474 extern PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 475 476 /*J 477 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 478 with an optional dynamic library name, for example 479 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 480 481 Level: beginner 482 483 .seealso: TSGLSetAcceptType(), TS 484 J*/ 485 #define TSGLAcceptType char* 486 #define TSGLACCEPT_ALWAYS "always" 487 488 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 489 extern PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 490 491 /*MC 492 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 493 494 Synopsis: 495 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 496 497 Not Collective 498 499 Input Parameters: 500 + name_scheme - name of user-defined acceptance scheme 501 . path - path (either absolute or relative) the library containing this scheme 502 . name_create - name of routine to create method context 503 - routine_create - routine to create method context 504 505 Notes: 506 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 507 508 If dynamic libraries are used, then the fourth input argument (routine_create) 509 is ignored. 510 511 Sample usage: 512 .vb 513 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 514 "MySchemeCreate",MySchemeCreate); 515 .ve 516 517 Then, your scheme can be chosen with the procedural interface via 518 $ TSGLSetAcceptType(ts,"my_scheme") 519 or at runtime via the option 520 $ -ts_gl_accept_type my_scheme 521 522 Level: advanced 523 524 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 525 and others of the form ${any_environmental_variable} occuring in pathname will be 526 replaced with appropriate values. 527 528 .keywords: TSGL, TSGLAcceptType, register 529 530 .seealso: TSGLRegisterAll() 531 M*/ 532 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 533 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 534 #else 535 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 536 #endif 537 538 /*J 539 TSGLType - family of time integration method within the General Linear class 540 541 Level: beginner 542 543 .seealso: TSGLSetType(), TSGLRegister() 544 J*/ 545 #define TSGLType char* 546 #define TSGL_IRKS "irks" 547 548 /*MC 549 TSGLRegisterDynamic - adds a TSGL implementation 550 551 Synopsis: 552 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 553 554 Not Collective 555 556 Input Parameters: 557 + name_scheme - name of user-defined general linear scheme 558 . path - path (either absolute or relative) the library containing this scheme 559 . name_create - name of routine to create method context 560 - routine_create - routine to create method context 561 562 Notes: 563 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 564 565 If dynamic libraries are used, then the fourth input argument (routine_create) 566 is ignored. 567 568 Sample usage: 569 .vb 570 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 571 "MySchemeCreate",MySchemeCreate); 572 .ve 573 574 Then, your scheme can be chosen with the procedural interface via 575 $ TSGLSetType(ts,"my_scheme") 576 or at runtime via the option 577 $ -ts_gl_type my_scheme 578 579 Level: advanced 580 581 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 582 and others of the form ${any_environmental_variable} occuring in pathname will be 583 replaced with appropriate values. 584 585 .keywords: TSGL, register 586 587 .seealso: TSGLRegisterAll() 588 M*/ 589 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 590 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 591 #else 592 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 593 #endif 594 595 extern PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 596 extern PetscErrorCode TSGLRegisterAll(const char[]); 597 extern PetscErrorCode TSGLRegisterDestroy(void); 598 extern PetscErrorCode TSGLInitializePackage(const char[]); 599 extern PetscErrorCode TSGLFinalizePackage(void); 600 extern PetscErrorCode TSGLSetType(TS,const TSGLType); 601 extern PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 602 extern PetscErrorCode TSGLSetAcceptType(TS,const TSGLAcceptType); 603 604 #define TSARKIMEXType char* 605 #define TSARKIMEX2D "2d" 606 #define TSARKIMEX2E "2e" 607 #define TSARKIMEXPRSSP2 "prssp2" 608 #define TSARKIMEX3 "3" 609 #define TSARKIMEXBPR3 "bpr3" 610 #define TSARKIMEXARS443 "ars443" 611 #define TSARKIMEX4 "4" 612 #define TSARKIMEX5 "5" 613 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*); 614 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType); 615 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 616 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]); 617 extern PetscErrorCode TSARKIMEXFinalizePackage(void); 618 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 619 extern PetscErrorCode TSARKIMEXRegisterDestroy(void); 620 extern PetscErrorCode TSARKIMEXRegisterAll(void); 621 622 #define TSRosWType char* 623 #define TSROSW2M "2m" 624 #define TSROSW2P "2p" 625 #define TSROSWRA3PW "ra3pw" 626 #define TSROSWRA34PW2 "ra34pw2" 627 #define TSROSWRODAS3 "rodas3" 628 #define TSROSWSANDU3 "sandu3" 629 #define TSROSWASSP3P3S1C "assp3p3s1c" 630 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 631 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 632 #define TSROSWARK3 "ark3" 633 634 extern PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*); 635 extern PetscErrorCode TSRosWSetType(TS ts,const TSRosWType); 636 extern PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 637 extern PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[]); 638 extern PetscErrorCode TSRosWFinalizePackage(void); 639 extern PetscErrorCode TSRosWInitializePackage(const char path[]); 640 extern PetscErrorCode TSRosWRegisterDestroy(void); 641 extern PetscErrorCode TSRosWRegisterAll(void); 642 643 /* 644 PETSc interface to Sundials 645 */ 646 #ifdef PETSC_HAVE_SUNDIALS 647 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 648 extern const char *TSSundialsLmmTypes[]; 649 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 650 extern const char *TSSundialsGramSchmidtTypes[]; 651 extern PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 652 extern PetscErrorCode TSSundialsGetPC(TS,PC*); 653 extern PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 654 extern PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 655 extern PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 656 extern PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 657 extern PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 658 extern PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 659 extern PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 660 extern PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 661 extern PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 662 extern PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 663 #endif 664 665 extern PetscErrorCode TSRKSetTolerance(TS,PetscReal); 666 667 extern PetscErrorCode TSThetaSetTheta(TS,PetscReal); 668 extern PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 669 extern PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 670 extern PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 671 672 extern PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 673 extern PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 674 extern PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 675 extern PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 676 extern PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 677 678 extern PetscErrorCode TSSetDM(TS,DM); 679 extern PetscErrorCode TSGetDM(TS,DM*); 680 681 extern PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 682 extern PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 683 684 PETSC_EXTERN_CXX_END 685 #endif 686