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