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