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