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