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