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