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