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