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