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