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