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