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 #define TSEIMEX "eimex" 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 TSEquationType - type of TS problem that is solved 55 56 Level: beginner 57 58 Developer Notes: this must match finclude/petscts.h 59 60 Supported types are: 61 TS_EQ_UNSPECIFIED (default) 62 TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0 63 TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0 64 65 .seealso: TSGetEquationType(), TSSetEquationType() 66 E*/ 67 typedef enum { 68 TS_EQ_UNSPECIFIED = -1, 69 TS_EQ_EXPLICIT = 0, 70 TS_EQ_ODE_EXPLICIT = 1, 71 TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 72 TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 73 TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 74 TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 75 TS_EQ_IMPLICIT = 1000, 76 TS_EQ_ODE_IMPLICIT = 1001, 77 TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 78 TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 79 TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 80 TS_EQ_DAE_IMPLICIT_INDEXHI = 1500 81 } TSEquationType; 82 PETSC_EXTERN const char *const*TSEquationTypes; 83 84 /*E 85 TSConvergedReason - reason a TS method has converged or not 86 87 Level: beginner 88 89 Developer Notes: this must match finclude/petscts.h 90 91 Each reason has its own manual page. 92 93 .seealso: TSGetConvergedReason() 94 E*/ 95 typedef enum { 96 TS_CONVERGED_ITERATING = 0, 97 TS_CONVERGED_TIME = 1, 98 TS_CONVERGED_ITS = 2, 99 TS_CONVERGED_USER = 3, 100 TS_DIVERGED_NONLINEAR_SOLVE = -1, 101 TS_DIVERGED_STEP_REJECTED = -2 102 } TSConvergedReason; 103 PETSC_EXTERN const char *const*TSConvergedReasons; 104 /*MC 105 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 106 107 Level: beginner 108 109 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt() 110 M*/ 111 112 /*MC 113 TS_CONVERGED_TIME - the final time was reached 114 115 Level: beginner 116 117 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime() 118 M*/ 119 120 /*MC 121 TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time 122 123 Level: beginner 124 125 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration() 126 M*/ 127 /*MC 128 TS_CONVERGED_USER - user requested termination 129 130 Level: beginner 131 132 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 133 M*/ 134 135 /*MC 136 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 137 138 Level: beginner 139 140 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason() 141 M*/ 142 143 /*MC 144 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 145 146 Level: beginner 147 148 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt() 149 M*/ 150 151 /*E 152 TSExactFinalTimeOption - option for handling of final time step 153 154 Level: beginner 155 156 Developer Notes: this must match finclude/petscts.h 157 158 $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded 159 $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time 160 $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time 161 .seealso: TSGetConvergedReason(), TSSetExactFinalTime() 162 163 E*/ 164 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption; 165 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 166 167 168 /* Logging support */ 169 PETSC_EXTERN PetscClassId TS_CLASSID; 170 PETSC_EXTERN PetscClassId DMTS_CLASSID; 171 172 PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 173 174 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 175 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 176 177 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 178 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 179 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 180 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 181 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 182 183 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 184 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 185 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 186 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 187 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 188 PETSC_EXTERN PetscErrorCode TSReset(TS); 189 190 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 191 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 192 193 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 194 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 195 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 196 197 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 198 199 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 200 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 201 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 202 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 203 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*); 204 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 205 206 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 207 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 208 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 209 210 PETSC_EXTERN PetscErrorCode TSStep(TS); 211 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 212 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 213 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*); 214 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType); 215 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 216 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 217 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 218 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 219 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 220 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 221 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 222 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 223 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 224 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 225 226 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 227 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 228 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 229 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 230 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 231 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 232 233 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 234 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 235 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 236 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 237 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 238 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 239 240 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 241 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 242 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*),void*); 243 244 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 245 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 246 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 247 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 248 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 249 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 250 251 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 252 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 253 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 254 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 255 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 256 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec); 257 258 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 259 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 260 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 261 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 262 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 263 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 264 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 265 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 266 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 267 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 268 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 269 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 270 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 271 272 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 273 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*); 274 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 275 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 276 277 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 278 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *); 279 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 280 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 281 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 282 283 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 284 285 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 286 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 287 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 288 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 289 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 290 291 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 292 293 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 294 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 295 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 296 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 297 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 298 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 299 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 300 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 301 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 302 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 303 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*),void*); 304 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**TSForcingFunction)(TS,PetscReal,Vec,void*),void**); 305 306 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 307 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 308 309 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 310 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*); 311 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 312 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*); 313 314 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 315 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *); 316 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 317 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *); 318 319 typedef struct { 320 Vec ray; 321 VecScatter scatter; 322 PetscViewer viewer; 323 } TSMonitorDMDARayCtx; 324 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**); 325 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*); 326 327 328 /* Dynamic creation and loading functions */ 329 PETSC_EXTERN PetscFunctionList TSList; 330 PETSC_EXTERN PetscBool TSRegisterAllCalled; 331 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 332 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 333 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 334 PETSC_EXTERN PetscErrorCode TSRegisterAll(void); 335 336 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 337 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 338 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 339 340 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 341 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 342 343 #define TS_FILE_CLASSID 1211225 344 345 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 346 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 347 348 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 349 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 350 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 351 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 352 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 353 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 354 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 355 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 356 357 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 358 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 359 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 360 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 361 362 /*J 363 TSSSPType - string with the name of TSSSP scheme. 364 365 Level: beginner 366 367 .seealso: TSSSPSetType(), TS 368 J*/ 369 typedef const char* TSSSPType; 370 #define TSSSPRKS2 "rks2" 371 #define TSSSPRKS3 "rks3" 372 #define TSSSPRK104 "rk104" 373 374 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 375 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 376 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 377 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 378 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 379 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 380 381 /*S 382 TSAdapt - Abstract object that manages time-step adaptivity 383 384 Level: beginner 385 386 .seealso: TS, TSAdaptCreate(), TSAdaptType 387 S*/ 388 typedef struct _p_TSAdapt *TSAdapt; 389 390 /*E 391 TSAdaptType - String with the name of TSAdapt scheme or the creation function 392 with an optional dynamic library name, for example 393 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 394 395 Level: beginner 396 397 .seealso: TSAdaptSetType(), TS 398 E*/ 399 typedef const char *TSAdaptType; 400 #define TSADAPTBASIC "basic" 401 #define TSADAPTNONE "none" 402 #define TSADAPTCFL "cfl" 403 404 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*); 405 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt)); 406 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(void); 407 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 408 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 409 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 410 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 411 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 412 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 413 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 414 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 415 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 416 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 417 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 418 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 419 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 420 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 421 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 422 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 423 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 424 425 /*S 426 TSGLAdapt - Abstract object that manages time-step adaptivity 427 428 Level: beginner 429 430 Developer Notes: 431 This functionality should be replaced by the TSAdapt. 432 433 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 434 S*/ 435 typedef struct _p_TSGLAdapt *TSGLAdapt; 436 437 /*J 438 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 439 with an optional dynamic library name, for example 440 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 441 442 Level: beginner 443 444 .seealso: TSGLAdaptSetType(), TS 445 J*/ 446 typedef const char *TSGLAdaptType; 447 #define TSGLADAPT_NONE "none" 448 #define TSGLADAPT_SIZE "size" 449 #define TSGLADAPT_BOTH "both" 450 451 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt)); 452 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(void); 453 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void); 454 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 455 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 456 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 457 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 458 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 459 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 460 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 461 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 462 463 /*J 464 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 465 with an optional dynamic library name, for example 466 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 467 468 Level: beginner 469 470 .seealso: TSGLSetAcceptType(), TS 471 J*/ 472 typedef const char *TSGLAcceptType; 473 #define TSGLACCEPT_ALWAYS "always" 474 475 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 476 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction); 477 478 /*J 479 TSGLType - family of time integration method within the General Linear class 480 481 Level: beginner 482 483 .seealso: TSGLSetType(), TSGLRegister() 484 J*/ 485 typedef const char* TSGLType; 486 #define TSGL_IRKS "irks" 487 488 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS)); 489 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(void); 490 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void); 491 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 492 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 493 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 494 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 495 496 /*J 497 TSEIMEXType - String with the name of an Extrapolated IMEX method. 498 499 Level: beginner 500 501 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister() 502 J*/ 503 #define TSEIMEXType char* 504 505 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt); 506 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt); 507 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool); 508 509 /*J 510 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 511 512 Level: beginner 513 514 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 515 J*/ 516 typedef const char* TSARKIMEXType; 517 #define TSARKIMEX1BEE "1bee" 518 #define TSARKIMEXA2 "a2" 519 #define TSARKIMEXL2 "l2" 520 #define TSARKIMEXARS122 "ars122" 521 #define TSARKIMEX2C "2c" 522 #define TSARKIMEX2D "2d" 523 #define TSARKIMEX2E "2e" 524 #define TSARKIMEXPRSSP2 "prssp2" 525 #define TSARKIMEX3 "3" 526 #define TSARKIMEXBPR3 "bpr3" 527 #define TSARKIMEXARS443 "ars443" 528 #define TSARKIMEX4 "4" 529 #define TSARKIMEX5 "5" 530 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 531 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 532 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 533 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[]); 534 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 535 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 536 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 537 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void); 538 539 /*J 540 TSRosWType - String with the name of a Rosenbrock-W method. 541 542 Level: beginner 543 544 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 545 J*/ 546 typedef const char* TSRosWType; 547 #define TSROSW2M "2m" 548 #define TSROSW2P "2p" 549 #define TSROSWRA3PW "ra3pw" 550 #define TSROSWRA34PW2 "ra34pw2" 551 #define TSROSWRODAS3 "rodas3" 552 #define TSROSWSANDU3 "sandu3" 553 #define TSROSWASSP3P3S1C "assp3p3s1c" 554 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 555 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 556 #define TSROSWARK3 "ark3" 557 #define TSROSWTHETA1 "theta1" 558 #define TSROSWTHETA2 "theta2" 559 #define TSROSWGRK4T "grk4t" 560 #define TSROSWSHAMP4 "shamp4" 561 #define TSROSWVELDD4 "veldd4" 562 #define TSROSW4L "4l" 563 564 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 565 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 566 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 567 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 568 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 569 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 570 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 571 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 572 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void); 573 574 /* 575 PETSc interface to Sundials 576 */ 577 #ifdef PETSC_HAVE_SUNDIALS 578 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 579 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 580 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 581 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 582 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 583 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 584 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 585 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 586 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 587 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 588 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 589 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 590 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 591 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 592 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 593 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 594 #endif 595 596 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal); 597 598 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 599 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 600 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 601 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 602 603 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 604 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 605 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 606 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 607 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 608 609 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 610 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 611 612 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 613 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 614 615 #endif 616