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