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_CONVERGED_USER = 3, 69 TS_DIVERGED_NONLINEAR_SOLVE = -1, 70 TS_DIVERGED_STEP_REJECTED = -2 71 } TSConvergedReason; 72 PETSC_EXTERN const char *const*TSConvergedReasons; 73 /*MC 74 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 75 76 Level: beginner 77 78 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt() 79 M*/ 80 81 /*MC 82 TS_CONVERGED_TIME - the final time was reached 83 84 Level: beginner 85 86 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt(), TSSetDuration(), TSGetSolveTime() 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(), TSGetConvergedReason(), TSGetTSAdapt(), TSSetDuration() 95 M*/ 96 /*MC 97 TS_CONVERGED_USER - user requested termination 98 99 Level: beginner 100 101 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 102 M*/ 103 104 /*MC 105 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 106 107 Level: beginner 108 109 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt(), TSGetSNES(), SNESGetConvergedReason() 110 M*/ 111 112 /*MC 113 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 114 115 Level: beginner 116 117 .seealso: TSSolve(), TSGetConvergedReason(), TSGetTSAdapt() 118 M*/ 119 120 /*E 121 TSExactFinalTimeOption - option for handling of final time step 122 123 Level: beginner 124 125 Developer Notes: this must match finclude/petscts.h 126 127 $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded 128 $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time 129 $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time 130 .seealso: TSGetConvergedReason() 131 E*/ 132 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption; 133 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 134 135 136 /* Logging support */ 137 PETSC_EXTERN PetscClassId TS_CLASSID; 138 PETSC_EXTERN PetscClassId DMTS_CLASSID; 139 140 PETSC_EXTERN PetscErrorCode TSInitializePackage(const char[]); 141 142 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 143 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 144 145 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 146 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 147 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 148 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 149 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 150 151 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 152 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 153 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 154 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 155 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 156 PETSC_EXTERN PetscErrorCode TSReset(TS); 157 158 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 159 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 160 161 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 162 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 163 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 164 165 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 166 167 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 168 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 169 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 170 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 171 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 172 173 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 174 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 175 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 176 177 PETSC_EXTERN PetscErrorCode TSStep(TS); 178 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 179 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 180 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 181 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 182 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 183 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 184 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 185 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 186 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 187 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 188 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 189 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 190 191 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 192 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 193 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 194 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 195 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 196 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 197 198 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 199 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 200 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 201 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 202 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 203 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 204 205 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 206 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 207 208 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 209 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 210 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 211 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 212 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 213 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 214 215 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 216 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 217 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 218 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 219 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 220 221 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 222 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 223 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 224 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 225 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 226 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 227 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 228 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 229 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 230 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 231 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 232 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 233 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 234 235 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 236 PETSC_EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 237 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 238 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 239 240 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 241 PETSC_EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 242 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 243 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 244 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 245 246 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 247 248 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 249 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 250 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 251 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 252 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 253 254 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 255 256 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 257 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 258 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 259 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 260 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 261 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 262 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 263 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 264 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 265 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 266 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 267 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 268 269 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 270 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*); 271 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 272 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*); 273 274 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 275 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *); 276 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 277 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *); 278 279 /* Dynamic creation and loading functions */ 280 PETSC_EXTERN PetscFList TSList; 281 PETSC_EXTERN PetscBool TSRegisterAllCalled; 282 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 283 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 284 PETSC_EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 285 PETSC_EXTERN PetscErrorCode TSRegisterAll(const char[]); 286 PETSC_EXTERN PetscErrorCode TSRegisterDestroy(void); 287 288 /*MC 289 TSRegisterDynamic - Adds a creation method to the TS package. 290 291 Synopsis: 292 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 293 294 Not Collective 295 296 Input Parameters: 297 + name - The name of a new user-defined creation routine 298 . path - The path (either absolute or relative) of the library containing this routine 299 . func_name - The name of the creation routine 300 - create_func - The creation routine itself 301 302 Notes: 303 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 304 305 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 306 307 Sample usage: 308 .vb 309 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 310 .ve 311 312 Then, your ts type can be chosen with the procedural interface via 313 .vb 314 TS ts; 315 TSCreate(MPI_Comm, &ts); 316 TSSetType(ts, "my_ts") 317 .ve 318 or at runtime via the option 319 .vb 320 -ts_type my_ts 321 .ve 322 323 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 324 If your function is not being put into a shared library then use TSRegister() instead 325 326 Level: advanced 327 328 .keywords: TS, register 329 .seealso: TSRegisterAll(), TSRegisterDestroy() 330 M*/ 331 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 332 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 333 #else 334 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 335 #endif 336 337 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 338 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 339 340 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 341 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 342 343 #define TS_FILE_CLASSID 1211225 344 345 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 346 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 347 348 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 349 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 350 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 351 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 352 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 353 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 354 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 355 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 356 357 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 358 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 359 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 360 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 361 362 /*J 363 TSSSPType - string with the name of TSSSP scheme. 364 365 Level: beginner 366 367 .seealso: TSSSPSetType(), TS 368 J*/ 369 typedef const char* TSSSPType; 370 #define TSSSPRKS2 "rks2" 371 #define TSSSPRKS3 "rks3" 372 #define TSSSPRK104 "rk104" 373 374 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 375 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 376 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 377 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 378 379 /*S 380 TSAdapt - Abstract object that manages time-step adaptivity 381 382 Level: beginner 383 384 .seealso: TS, TSAdaptCreate(), TSAdaptType 385 S*/ 386 typedef struct _p_TSAdapt *TSAdapt; 387 388 /*E 389 TSAdaptType - String with the name of TSAdapt scheme or the creation function 390 with an optional dynamic library name, for example 391 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 392 393 Level: beginner 394 395 .seealso: TSAdaptSetType(), TS 396 E*/ 397 typedef const char *TSAdaptType; 398 #define TSADAPTBASIC "basic" 399 #define TSADAPTNONE "none" 400 #define TSADAPTCFL "cfl" 401 402 /*MC 403 TSAdaptRegisterDynamic - adds a TSAdapt implementation 404 405 Synopsis: 406 PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 407 408 Not Collective 409 410 Input Parameters: 411 + name_scheme - name of user-defined adaptivity scheme 412 . path - path (either absolute or relative) the library containing this scheme 413 . name_create - name of routine to create method context 414 - routine_create - routine to create method context 415 416 Notes: 417 TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 418 419 If dynamic libraries are used, then the fourth input argument (routine_create) 420 is ignored. 421 422 Sample usage: 423 .vb 424 TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 425 "MySchemeCreate",MySchemeCreate); 426 .ve 427 428 Then, your scheme can be chosen with the procedural interface via 429 $ TSAdaptSetType(ts,"my_scheme") 430 or at runtime via the option 431 $ -ts_adapt_type my_scheme 432 433 Level: advanced 434 435 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 436 and others of the form ${any_environmental_variable} occuring in pathname will be 437 replaced with appropriate values. 438 439 .keywords: TSAdapt, register 440 441 .seealso: TSAdaptRegisterAll() 442 M*/ 443 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 444 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0) 445 #else 446 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d) 447 #endif 448 449 PETSC_EXTERN PetscErrorCode TSGetTSAdapt(TS,TSAdapt*); 450 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt)); 451 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]); 452 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void); 453 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]); 454 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 455 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 456 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 457 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 458 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 459 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 460 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 461 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 462 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 463 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 464 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 465 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 466 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 467 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 468 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 469 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 470 471 /*S 472 TSGLAdapt - Abstract object that manages time-step adaptivity 473 474 Level: beginner 475 476 Developer Notes: 477 This functionality should be replaced by the TSAdapt. 478 479 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 480 S*/ 481 typedef struct _p_TSGLAdapt *TSGLAdapt; 482 483 /*J 484 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 485 with an optional dynamic library name, for example 486 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 487 488 Level: beginner 489 490 .seealso: TSGLAdaptSetType(), TS 491 J*/ 492 typedef const char *TSGLAdaptType; 493 #define TSGLADAPT_NONE "none" 494 #define TSGLADAPT_SIZE "size" 495 #define TSGLADAPT_BOTH "both" 496 497 /*MC 498 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 499 500 Synopsis: 501 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 502 503 Not Collective 504 505 Input Parameters: 506 + name_scheme - name of user-defined adaptivity scheme 507 . path - path (either absolute or relative) the library containing this scheme 508 . name_create - name of routine to create method context 509 - routine_create - routine to create method context 510 511 Notes: 512 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 513 514 If dynamic libraries are used, then the fourth input argument (routine_create) 515 is ignored. 516 517 Sample usage: 518 .vb 519 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 520 "MySchemeCreate",MySchemeCreate); 521 .ve 522 523 Then, your scheme can be chosen with the procedural interface via 524 $ TSGLAdaptSetType(ts,"my_scheme") 525 or at runtime via the option 526 $ -ts_adapt_type my_scheme 527 528 Level: advanced 529 530 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 531 and others of the form ${any_environmental_variable} occuring in pathname will be 532 replaced with appropriate values. 533 534 .keywords: TSGLAdapt, register 535 536 .seealso: TSGLAdaptRegisterAll() 537 M*/ 538 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 539 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 540 #else 541 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 542 #endif 543 544 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 545 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]); 546 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void); 547 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]); 548 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 549 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 550 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 551 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 552 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 553 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 554 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 555 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 556 557 /*J 558 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 559 with an optional dynamic library name, for example 560 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 561 562 Level: beginner 563 564 .seealso: TSGLSetAcceptType(), TS 565 J*/ 566 typedef const char *TSGLAcceptType; 567 #define TSGLACCEPT_ALWAYS "always" 568 569 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 570 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 571 572 /*MC 573 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 574 575 Synopsis: 576 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 577 578 Not Collective 579 580 Input Parameters: 581 + name_scheme - name of user-defined acceptance scheme 582 . path - path (either absolute or relative) the library containing this scheme 583 . name_create - name of routine to create method context 584 - routine_create - routine to create method context 585 586 Notes: 587 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 588 589 If dynamic libraries are used, then the fourth input argument (routine_create) 590 is ignored. 591 592 Sample usage: 593 .vb 594 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 595 "MySchemeCreate",MySchemeCreate); 596 .ve 597 598 Then, your scheme can be chosen with the procedural interface via 599 $ TSGLSetAcceptType(ts,"my_scheme") 600 or at runtime via the option 601 $ -ts_gl_accept_type my_scheme 602 603 Level: advanced 604 605 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 606 and others of the form ${any_environmental_variable} occuring in pathname will be 607 replaced with appropriate values. 608 609 .keywords: TSGL, TSGLAcceptType, register 610 611 .seealso: TSGLRegisterAll() 612 M*/ 613 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 614 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 615 #else 616 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 617 #endif 618 619 /*J 620 TSGLType - family of time integration method within the General Linear class 621 622 Level: beginner 623 624 .seealso: TSGLSetType(), TSGLRegister() 625 J*/ 626 typedef const char* TSGLType; 627 #define TSGL_IRKS "irks" 628 629 /*MC 630 TSGLRegisterDynamic - adds a TSGL implementation 631 632 Synopsis: 633 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 634 635 Not Collective 636 637 Input Parameters: 638 + name_scheme - name of user-defined general linear scheme 639 . path - path (either absolute or relative) the library containing this scheme 640 . name_create - name of routine to create method context 641 - routine_create - routine to create method context 642 643 Notes: 644 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 645 646 If dynamic libraries are used, then the fourth input argument (routine_create) 647 is ignored. 648 649 Sample usage: 650 .vb 651 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 652 "MySchemeCreate",MySchemeCreate); 653 .ve 654 655 Then, your scheme can be chosen with the procedural interface via 656 $ TSGLSetType(ts,"my_scheme") 657 or at runtime via the option 658 $ -ts_gl_type my_scheme 659 660 Level: advanced 661 662 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 663 and others of the form ${any_environmental_variable} occuring in pathname will be 664 replaced with appropriate values. 665 666 .keywords: TSGL, register 667 668 .seealso: TSGLRegisterAll() 669 M*/ 670 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 671 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 672 #else 673 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 674 #endif 675 676 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 677 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]); 678 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void); 679 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]); 680 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 681 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 682 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 683 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 684 685 /*J 686 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 687 688 Level: beginner 689 690 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 691 J*/ 692 typedef const char* TSARKIMEXType; 693 #define TSARKIMEXA2 "a2" 694 #define TSARKIMEXL2 "l2" 695 #define TSARKIMEXARS122 "ars122" 696 #define TSARKIMEX2C "2c" 697 #define TSARKIMEX2D "2d" 698 #define TSARKIMEX2E "2e" 699 #define TSARKIMEXPRSSP2 "prssp2" 700 #define TSARKIMEX3 "3" 701 #define TSARKIMEXBPR3 "bpr3" 702 #define TSARKIMEXARS443 "ars443" 703 #define TSARKIMEX4 "4" 704 #define TSARKIMEX5 "5" 705 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 706 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 707 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 708 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[]); 709 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 710 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 711 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 712 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void); 713 714 /*J 715 TSRosWType - String with the name of a Rosenbrock-W method. 716 717 Level: beginner 718 719 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 720 J*/ 721 typedef const char* TSRosWType; 722 #define TSROSW2M "2m" 723 #define TSROSW2P "2p" 724 #define TSROSWRA3PW "ra3pw" 725 #define TSROSWRA34PW2 "ra34pw2" 726 #define TSROSWRODAS3 "rodas3" 727 #define TSROSWSANDU3 "sandu3" 728 #define TSROSWASSP3P3S1C "assp3p3s1c" 729 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 730 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 731 #define TSROSWARK3 "ark3" 732 #define TSROSWTHETA1 "theta1" 733 #define TSROSWTHETA2 "theta2" 734 #define TSROSWGRK4T "grk4t" 735 #define TSROSWSHAMP4 "shamp4" 736 #define TSROSWVELDD4 "veldd4" 737 #define TSROSW4L "4l" 738 739 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 740 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 741 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 742 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 743 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 744 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 745 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]); 746 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 747 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void); 748 749 /* 750 PETSc interface to Sundials 751 */ 752 #ifdef PETSC_HAVE_SUNDIALS 753 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 754 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 755 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 756 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 757 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 758 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 759 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 760 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 761 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 762 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 763 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 764 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 765 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 766 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 767 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 768 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 769 #endif 770 771 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal); 772 773 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 774 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 775 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 776 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 777 778 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 779 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 780 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 781 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 782 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 783 784 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 785 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 786 787 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 788 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 789 790 #endif 791