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