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 extern PetscErrorCode TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 183 extern PetscErrorCode TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 184 185 extern PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 186 extern PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 187 extern PetscErrorCode TSPreStep(TS); 188 extern PetscErrorCode TSPostStep(TS); 189 extern PetscErrorCode TSSetRetainStages(TS,PetscBool); 190 extern PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 191 extern PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 192 extern PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 193 extern PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 194 extern PetscErrorCode TSGetCFLTime(TS,PetscReal*); 195 196 extern PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 197 extern PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 198 extern PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 199 extern PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 200 201 extern PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 202 extern PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 203 extern PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 204 extern PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 205 extern PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 206 207 extern PetscErrorCode TSPythonSetType(TS,const char[]); 208 209 extern PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 210 extern PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 211 extern PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 212 extern PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 213 214 extern PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 215 216 /* Dynamic creation and loading functions */ 217 extern PetscFList TSList; 218 extern PetscBool TSRegisterAllCalled; 219 extern PetscErrorCode TSGetType(TS,const TSType*); 220 extern PetscErrorCode TSSetType(TS,const TSType); 221 extern PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 222 extern PetscErrorCode TSRegisterAll(const char[]); 223 extern PetscErrorCode TSRegisterDestroy(void); 224 225 /*MC 226 TSRegisterDynamic - Adds a creation method to the TS package. 227 228 Synopsis: 229 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 230 231 Not Collective 232 233 Input Parameters: 234 + name - The name of a new user-defined creation routine 235 . path - The path (either absolute or relative) of the library containing this routine 236 . func_name - The name of the creation routine 237 - create_func - The creation routine itself 238 239 Notes: 240 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 241 242 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 243 244 Sample usage: 245 .vb 246 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 247 .ve 248 249 Then, your ts type can be chosen with the procedural interface via 250 .vb 251 TS ts; 252 TSCreate(MPI_Comm, &ts); 253 TSSetType(ts, "my_ts") 254 .ve 255 or at runtime via the option 256 .vb 257 -ts_type my_ts 258 .ve 259 260 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 261 If your function is not being put into a shared library then use TSRegister() instead 262 263 Level: advanced 264 265 .keywords: TS, register 266 .seealso: TSRegisterAll(), TSRegisterDestroy() 267 M*/ 268 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 269 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 270 #else 271 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 272 #endif 273 274 extern PetscErrorCode TSGetSNES(TS,SNES*); 275 extern PetscErrorCode TSGetKSP(TS,KSP*); 276 277 extern PetscErrorCode TSView(TS,PetscViewer); 278 279 extern PetscErrorCode TSSetApplicationContext(TS,void *); 280 extern PetscErrorCode TSGetApplicationContext(TS,void *); 281 282 extern PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *); 283 extern PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *); 284 extern PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*); 285 286 /*J 287 TSSSPType - string with the name of TSSSP scheme. 288 289 Level: beginner 290 291 .seealso: TSSSPSetType(), TS 292 J*/ 293 #define TSSSPType char* 294 #define TSSSPRKS2 "rks2" 295 #define TSSSPRKS3 "rks3" 296 #define TSSSPRK104 "rk104" 297 298 extern PetscErrorCode TSSSPSetType(TS,const TSSSPType); 299 extern PetscErrorCode TSSSPGetType(TS,const TSSSPType*); 300 extern PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 301 extern PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 302 303 /*S 304 TSAdapt - Abstract object that manages time-step adaptivity 305 306 Level: beginner 307 308 .seealso: TS, TSAdaptCreate(), TSAdaptType 309 S*/ 310 typedef struct _p_TSAdapt *TSAdapt; 311 312 /*E 313 TSAdaptType - String with the name of TSAdapt scheme or the creation function 314 with an optional dynamic library name, for example 315 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 316 317 Level: beginner 318 319 .seealso: TSAdaptSetType(), TS 320 E*/ 321 #define TSAdaptType char* 322 #define TSADAPTBASIC "basic" 323 #define TSADAPTNONE "none" 324 #define TSADAPTCFL "cfl" 325 326 /*MC 327 TSAdaptRegisterDynamic - adds a TSAdapt implementation 328 329 Synopsis: 330 PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 331 332 Not Collective 333 334 Input Parameters: 335 + name_scheme - name of user-defined adaptivity scheme 336 . path - path (either absolute or relative) the library containing this scheme 337 . name_create - name of routine to create method context 338 - routine_create - routine to create method context 339 340 Notes: 341 TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 342 343 If dynamic libraries are used, then the fourth input argument (routine_create) 344 is ignored. 345 346 Sample usage: 347 .vb 348 TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 349 "MySchemeCreate",MySchemeCreate); 350 .ve 351 352 Then, your scheme can be chosen with the procedural interface via 353 $ TSAdaptSetType(ts,"my_scheme") 354 or at runtime via the option 355 $ -ts_adapt_type my_scheme 356 357 Level: advanced 358 359 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 360 and others of the form ${any_environmental_variable} occuring in pathname will be 361 replaced with appropriate values. 362 363 .keywords: TSAdapt, register 364 365 .seealso: TSAdaptRegisterAll() 366 M*/ 367 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 368 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0) 369 #else 370 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d) 371 #endif 372 373 extern PetscErrorCode TSGetAdapt(TS,TSAdapt*); 374 extern PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt)); 375 extern PetscErrorCode TSAdaptRegisterAll(const char[]); 376 extern PetscErrorCode TSAdaptRegisterDestroy(void); 377 extern PetscErrorCode TSAdaptInitializePackage(const char[]); 378 extern PetscErrorCode TSAdaptFinalizePackage(void); 379 extern PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 380 extern PetscErrorCode TSAdaptSetType(TSAdapt,const TSAdaptType); 381 extern PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 382 extern PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 383 extern PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 384 extern PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 385 extern PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 386 extern PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 387 extern PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 388 extern PetscErrorCode TSAdaptDestroy(TSAdapt*); 389 extern PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 390 extern PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 391 392 /*S 393 TSGLAdapt - Abstract object that manages time-step adaptivity 394 395 Level: beginner 396 397 Developer Notes: 398 This functionality should be replaced by the TSAdapt. 399 400 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 401 S*/ 402 typedef struct _p_TSGLAdapt *TSGLAdapt; 403 404 /*J 405 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 406 with an optional dynamic library name, for example 407 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 408 409 Level: beginner 410 411 .seealso: TSGLAdaptSetType(), TS 412 J*/ 413 #define TSGLAdaptType char* 414 #define TSGLADAPT_NONE "none" 415 #define TSGLADAPT_SIZE "size" 416 #define TSGLADAPT_BOTH "both" 417 418 /*MC 419 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 420 421 Synopsis: 422 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 423 424 Not Collective 425 426 Input Parameters: 427 + name_scheme - name of user-defined adaptivity scheme 428 . path - path (either absolute or relative) the library containing this scheme 429 . name_create - name of routine to create method context 430 - routine_create - routine to create method context 431 432 Notes: 433 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 434 435 If dynamic libraries are used, then the fourth input argument (routine_create) 436 is ignored. 437 438 Sample usage: 439 .vb 440 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 441 "MySchemeCreate",MySchemeCreate); 442 .ve 443 444 Then, your scheme can be chosen with the procedural interface via 445 $ TSGLAdaptSetType(ts,"my_scheme") 446 or at runtime via the option 447 $ -ts_adapt_type my_scheme 448 449 Level: advanced 450 451 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 452 and others of the form ${any_environmental_variable} occuring in pathname will be 453 replaced with appropriate values. 454 455 .keywords: TSGLAdapt, register 456 457 .seealso: TSGLAdaptRegisterAll() 458 M*/ 459 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 460 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 461 #else 462 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 463 #endif 464 465 extern PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 466 extern PetscErrorCode TSGLAdaptRegisterAll(const char[]); 467 extern PetscErrorCode TSGLAdaptRegisterDestroy(void); 468 extern PetscErrorCode TSGLAdaptInitializePackage(const char[]); 469 extern PetscErrorCode TSGLAdaptFinalizePackage(void); 470 extern PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 471 extern PetscErrorCode TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType); 472 extern PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 473 extern PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 474 extern PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 475 extern PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 476 extern PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 477 478 /*J 479 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 480 with an optional dynamic library name, for example 481 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 482 483 Level: beginner 484 485 .seealso: TSGLSetAcceptType(), TS 486 J*/ 487 #define TSGLAcceptType char* 488 #define TSGLACCEPT_ALWAYS "always" 489 490 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 491 extern PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 492 493 /*MC 494 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 495 496 Synopsis: 497 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 498 499 Not Collective 500 501 Input Parameters: 502 + name_scheme - name of user-defined acceptance scheme 503 . path - path (either absolute or relative) the library containing this scheme 504 . name_create - name of routine to create method context 505 - routine_create - routine to create method context 506 507 Notes: 508 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 509 510 If dynamic libraries are used, then the fourth input argument (routine_create) 511 is ignored. 512 513 Sample usage: 514 .vb 515 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 516 "MySchemeCreate",MySchemeCreate); 517 .ve 518 519 Then, your scheme can be chosen with the procedural interface via 520 $ TSGLSetAcceptType(ts,"my_scheme") 521 or at runtime via the option 522 $ -ts_gl_accept_type my_scheme 523 524 Level: advanced 525 526 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 527 and others of the form ${any_environmental_variable} occuring in pathname will be 528 replaced with appropriate values. 529 530 .keywords: TSGL, TSGLAcceptType, register 531 532 .seealso: TSGLRegisterAll() 533 M*/ 534 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 535 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 536 #else 537 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 538 #endif 539 540 /*J 541 TSGLType - family of time integration method within the General Linear class 542 543 Level: beginner 544 545 .seealso: TSGLSetType(), TSGLRegister() 546 J*/ 547 #define TSGLType char* 548 #define TSGL_IRKS "irks" 549 550 /*MC 551 TSGLRegisterDynamic - adds a TSGL implementation 552 553 Synopsis: 554 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 555 556 Not Collective 557 558 Input Parameters: 559 + name_scheme - name of user-defined general linear scheme 560 . path - path (either absolute or relative) the library containing this scheme 561 . name_create - name of routine to create method context 562 - routine_create - routine to create method context 563 564 Notes: 565 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 566 567 If dynamic libraries are used, then the fourth input argument (routine_create) 568 is ignored. 569 570 Sample usage: 571 .vb 572 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 573 "MySchemeCreate",MySchemeCreate); 574 .ve 575 576 Then, your scheme can be chosen with the procedural interface via 577 $ TSGLSetType(ts,"my_scheme") 578 or at runtime via the option 579 $ -ts_gl_type my_scheme 580 581 Level: advanced 582 583 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 584 and others of the form ${any_environmental_variable} occuring in pathname will be 585 replaced with appropriate values. 586 587 .keywords: TSGL, register 588 589 .seealso: TSGLRegisterAll() 590 M*/ 591 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 592 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 593 #else 594 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 595 #endif 596 597 extern PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 598 extern PetscErrorCode TSGLRegisterAll(const char[]); 599 extern PetscErrorCode TSGLRegisterDestroy(void); 600 extern PetscErrorCode TSGLInitializePackage(const char[]); 601 extern PetscErrorCode TSGLFinalizePackage(void); 602 extern PetscErrorCode TSGLSetType(TS,const TSGLType); 603 extern PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 604 extern PetscErrorCode TSGLSetAcceptType(TS,const TSGLAcceptType); 605 606 #define TSARKIMEXType char* 607 #define TSARKIMEX2D "2d" 608 #define TSARKIMEX2E "2e" 609 #define TSARKIMEXPRSSP2 "prssp2" 610 #define TSARKIMEX3 "3" 611 #define TSARKIMEXBPR3 "bpr3" 612 #define TSARKIMEXARS443 "ars443" 613 #define TSARKIMEX4 "4" 614 #define TSARKIMEX5 "5" 615 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*); 616 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType); 617 extern PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 618 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[]); 619 extern PetscErrorCode TSARKIMEXFinalizePackage(void); 620 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 621 extern PetscErrorCode TSARKIMEXRegisterDestroy(void); 622 extern PetscErrorCode TSARKIMEXRegisterAll(void); 623 624 #define TSRosWType char* 625 #define TSROSW2M "2m" 626 #define TSROSW2P "2p" 627 #define TSROSWRA3PW "ra3pw" 628 #define TSROSWRA34PW2 "ra34pw2" 629 #define TSROSWRODAS3 "rodas3" 630 #define TSROSWSANDU3 "sandu3" 631 #define TSROSWASSP3P3S1C "assp3p3s1c" 632 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 633 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 634 #define TSROSWARK3 "ark3" 635 636 extern PetscErrorCode TSRosWGetType(TS ts,const TSRosWType*); 637 extern PetscErrorCode TSRosWSetType(TS ts,const TSRosWType); 638 extern PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 639 extern PetscErrorCode TSRosWRegister(const TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[]); 640 extern PetscErrorCode TSRosWFinalizePackage(void); 641 extern PetscErrorCode TSRosWInitializePackage(const char path[]); 642 extern PetscErrorCode TSRosWRegisterDestroy(void); 643 extern PetscErrorCode TSRosWRegisterAll(void); 644 645 /* 646 PETSc interface to Sundials 647 */ 648 #ifdef PETSC_HAVE_SUNDIALS 649 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 650 extern const char *TSSundialsLmmTypes[]; 651 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 652 extern const char *TSSundialsGramSchmidtTypes[]; 653 extern PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 654 extern PetscErrorCode TSSundialsGetPC(TS,PC*); 655 extern PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 656 extern PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 657 extern PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 658 extern PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 659 extern PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 660 extern PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 661 extern PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 662 extern PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 663 extern PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 664 extern PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 665 #endif 666 667 extern PetscErrorCode TSRKSetTolerance(TS,PetscReal); 668 669 extern PetscErrorCode TSThetaSetTheta(TS,PetscReal); 670 extern PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 671 extern PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 672 extern PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 673 674 extern PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 675 extern PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 676 extern PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 677 extern PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 678 extern PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 679 680 extern PetscErrorCode TSSetDM(TS,DM); 681 extern PetscErrorCode TSGetDM(TS,DM*); 682 683 extern PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 684 extern PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 685 686 PETSC_EXTERN_CXX_END 687 #endif 688