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