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