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 petsc-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 petsc-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 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 (time-steps) was reached prior to the final time 122 123 Level: beginner 124 125 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration() 126 M*/ 127 128 /*MC 129 TS_CONVERGED_USER - user requested termination 130 131 Level: beginner 132 133 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 134 M*/ 135 136 /*MC 137 TS_CONVERGED_EVENT - user requested termination on event detection 138 139 Level: beginner 140 141 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 142 M*/ 143 144 /*MC 145 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 146 147 Level: beginner 148 149 Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures. 150 151 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures() 152 M*/ 153 154 /*MC 155 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 156 157 Level: beginner 158 159 Notes: See TSSetMaxStepRejections() for how to allow more step rejections. 160 161 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections() 162 M*/ 163 164 /*E 165 TSExactFinalTimeOption - option for handling of final time step 166 167 Level: beginner 168 169 Developer Notes: this must match petsc-finclude/petscts.h 170 171 $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded 172 $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time 173 $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time 174 .seealso: TSGetConvergedReason(), TSSetExactFinalTime() 175 176 E*/ 177 typedef enum {TS_EXACTFINALTIME_STEPOVER=0,TS_EXACTFINALTIME_INTERPOLATE=1,TS_EXACTFINALTIME_MATCHSTEP=2} TSExactFinalTimeOption; 178 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 179 180 181 /* Logging support */ 182 PETSC_EXTERN PetscClassId TS_CLASSID; 183 PETSC_EXTERN PetscClassId DMTS_CLASSID; 184 185 PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 186 187 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 188 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 189 190 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 191 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 192 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 193 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 194 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 195 196 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 197 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 198 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 199 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 200 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 201 PETSC_EXTERN PetscErrorCode TSReset(TS); 202 203 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 204 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 205 206 PETSC_EXTERN PetscErrorCode TSAdjointSetSensitivity(TS,Vec*,PetscInt); 207 PETSC_EXTERN PetscErrorCode TSAdjointGetSensitivity(TS,Vec**,PetscInt*); 208 209 PETSC_EXTERN PetscErrorCode TSAdjointSetSensitivityP(TS,Vec*,PetscInt); 210 PETSC_EXTERN PetscErrorCode TSAdjointGetSensitivityP(TS,Vec**,PetscInt*); 211 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobianP(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*); 212 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobianP(TS,PetscReal,Vec,Mat); 213 PETSC_EXTERN PetscErrorCode TSAdjointSetCostIntegrand(TS,PetscInt,Vec,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),void*); 214 PETSC_EXTERN PetscErrorCode TSAdjointSetDRDYFunction(TS,Vec*,PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),void*); 215 PETSC_EXTERN PetscErrorCode TSAdjointSetDRDPFunction(TS,Vec*,PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),void*); 216 PETSC_EXTERN PetscErrorCode TSAdjointGetCostQuadrature(TS,Vec*); 217 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS,Vec); 218 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 219 220 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 221 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*); 222 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*); 223 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec); 224 225 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 226 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 227 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 228 229 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 230 231 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 232 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 233 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 234 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 235 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*); 236 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 237 238 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 239 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 240 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 241 242 PETSC_EXTERN PetscErrorCode TSStep(TS); 243 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 244 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 245 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*); 246 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType); 247 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 248 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 249 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 250 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 251 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 252 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 253 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 254 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 255 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 256 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 257 PETSC_EXTERN PetscErrorCode TSRollBack(TS); 258 259 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**); 260 PETSC_EXTERN PetscErrorCode TSSetCheckpoint(TS,PetscBool); 261 262 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 263 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 264 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 265 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 266 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 267 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 268 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*); 269 270 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 271 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*); 272 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 273 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 274 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 275 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 276 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool); 277 278 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 279 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 280 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*); 281 282 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 283 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 284 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 285 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 286 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 287 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 288 289 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 290 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*); 291 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 292 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 293 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 294 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec); 295 296 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 297 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 298 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*)); 299 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 300 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 301 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 302 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*); 303 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 304 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 305 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 306 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 307 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 308 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 309 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 310 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 311 312 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 313 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*); 314 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 315 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 316 317 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 318 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *); 319 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 320 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 321 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 322 323 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 324 325 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 326 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat); 327 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 328 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool); 329 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 330 331 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 332 333 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 334 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 335 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 336 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 337 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 338 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 339 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 340 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 341 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 342 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 343 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*); 344 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**); 345 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*); 346 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal); 347 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, void (**)(const PetscReal[], PetscScalar *, void *), void **); 348 349 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*); 350 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*); 351 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*); 352 353 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 354 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 355 356 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 357 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*); 358 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 359 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*); 360 361 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 362 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *); 363 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 364 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *); 365 366 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *); 367 368 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 369 typedef struct { 370 Vec ray; 371 VecScatter scatter; 372 PetscViewer viewer; 373 TSMonitorLGCtx lgctx; 374 } TSMonitorDMDARayCtx; 375 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**); 376 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*); 377 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*); 378 379 380 /* Dynamic creation and loading functions */ 381 PETSC_EXTERN PetscFunctionList TSList; 382 PETSC_EXTERN PetscBool TSRegisterAllCalled; 383 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 384 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 385 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 386 PETSC_EXTERN PetscErrorCode TSRegisterAll(void); 387 388 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 389 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 390 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 391 392 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 393 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 394 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 395 396 397 #define TS_FILE_CLASSID 1211225 398 399 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 400 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 401 402 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 403 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 404 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 405 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 406 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 407 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 408 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 409 410 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 411 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 412 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 413 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 414 415 PETSC_EXTERN PetscErrorCode TSSetEventMonitor(TS,PetscInt,PetscInt*,PetscBool*,PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar*,void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,void*),void*); 416 /*J 417 TSSSPType - string with the name of TSSSP scheme. 418 419 Level: beginner 420 421 .seealso: TSSSPSetType(), TS 422 J*/ 423 typedef const char* TSSSPType; 424 #define TSSSPRKS2 "rks2" 425 #define TSSSPRKS3 "rks3" 426 #define TSSSPRK104 "rk104" 427 428 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 429 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 430 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 431 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 432 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 433 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 434 435 /*S 436 TSAdapt - Abstract object that manages time-step adaptivity 437 438 Level: beginner 439 440 .seealso: TS, TSAdaptCreate(), TSAdaptType 441 S*/ 442 typedef struct _p_TSAdapt *TSAdapt; 443 444 /*E 445 TSAdaptType - String with the name of TSAdapt scheme. 446 447 Level: beginner 448 449 .seealso: TSAdaptSetType(), TS 450 E*/ 451 typedef const char *TSAdaptType; 452 #define TSADAPTBASIC "basic" 453 #define TSADAPTNONE "none" 454 #define TSADAPTCFL "cfl" 455 456 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*); 457 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt)); 458 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(void); 459 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 460 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 461 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 462 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 463 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 464 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 465 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 466 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 467 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 468 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 469 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 470 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 471 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptions*,TSAdapt); 472 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 473 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 474 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 475 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 476 477 /*S 478 TSGLAdapt - Abstract object that manages time-step adaptivity 479 480 Level: beginner 481 482 Developer Notes: 483 This functionality should be replaced by the TSAdapt. 484 485 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 486 S*/ 487 typedef struct _p_TSGLAdapt *TSGLAdapt; 488 489 /*J 490 TSGLAdaptType - String with the name of TSGLAdapt scheme 491 492 Level: beginner 493 494 .seealso: TSGLAdaptSetType(), TS 495 J*/ 496 typedef const char *TSGLAdaptType; 497 #define TSGLADAPT_NONE "none" 498 #define TSGLADAPT_SIZE "size" 499 #define TSGLADAPT_BOTH "both" 500 501 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt)); 502 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(void); 503 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void); 504 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 505 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 506 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 507 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 508 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 509 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 510 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptions*,TSGLAdapt); 511 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 512 513 /*J 514 TSGLAcceptType - String with the name of TSGLAccept scheme 515 516 Level: beginner 517 518 .seealso: TSGLSetAcceptType(), TS 519 J*/ 520 typedef const char *TSGLAcceptType; 521 #define TSGLACCEPT_ALWAYS "always" 522 523 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 524 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction); 525 526 /*J 527 TSGLType - family of time integration method within the General Linear class 528 529 Level: beginner 530 531 .seealso: TSGLSetType(), TSGLRegister() 532 J*/ 533 typedef const char* TSGLType; 534 #define TSGL_IRKS "irks" 535 536 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS)); 537 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(void); 538 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void); 539 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 540 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 541 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 542 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 543 544 /*J 545 TSEIMEXType - String with the name of an Extrapolated IMEX method. 546 547 Level: beginner 548 549 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister() 550 J*/ 551 #define TSEIMEXType char* 552 553 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt); 554 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt); 555 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool); 556 557 /*J 558 TSRKType - String with the name of a Runge-Kutta method. 559 560 Level: beginner 561 562 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister() 563 J*/ 564 typedef const char* TSRKType; 565 #define TSRK1FE "1fe" 566 #define TSRK2A "2a" 567 #define TSRK3 "3" 568 #define TSRK3BS "3bs" 569 #define TSRK4 "4" 570 #define TSRK5F "5f" 571 #define TSRK5DP "5dp" 572 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*); 573 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType); 574 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool); 575 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 576 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 577 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 578 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 579 PETSC_EXTERN PetscErrorCode TSRKRegisterAll(void); 580 581 /*J 582 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 583 584 Level: beginner 585 586 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 587 J*/ 588 typedef const char* TSARKIMEXType; 589 #define TSARKIMEX1BEE "1bee" 590 #define TSARKIMEXA2 "a2" 591 #define TSARKIMEXL2 "l2" 592 #define TSARKIMEXARS122 "ars122" 593 #define TSARKIMEX2C "2c" 594 #define TSARKIMEX2D "2d" 595 #define TSARKIMEX2E "2e" 596 #define TSARKIMEXPRSSP2 "prssp2" 597 #define TSARKIMEX3 "3" 598 #define TSARKIMEXBPR3 "bpr3" 599 #define TSARKIMEXARS443 "ars443" 600 #define TSARKIMEX4 "4" 601 #define TSARKIMEX5 "5" 602 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 603 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 604 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 605 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[]); 606 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 607 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 608 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 609 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void); 610 611 /*J 612 TSRosWType - String with the name of a Rosenbrock-W method. 613 614 Level: beginner 615 616 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 617 J*/ 618 typedef const char* TSRosWType; 619 #define TSROSW2M "2m" 620 #define TSROSW2P "2p" 621 #define TSROSWRA3PW "ra3pw" 622 #define TSROSWRA34PW2 "ra34pw2" 623 #define TSROSWRODAS3 "rodas3" 624 #define TSROSWSANDU3 "sandu3" 625 #define TSROSWASSP3P3S1C "assp3p3s1c" 626 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 627 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 628 #define TSROSWARK3 "ark3" 629 #define TSROSWTHETA1 "theta1" 630 #define TSROSWTHETA2 "theta2" 631 #define TSROSWGRK4T "grk4t" 632 #define TSROSWSHAMP4 "shamp4" 633 #define TSROSWVELDD4 "veldd4" 634 #define TSROSW4L "4l" 635 636 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 637 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 638 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 639 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 640 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 641 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 642 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 643 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 644 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void); 645 646 /* 647 PETSc interface to Sundials 648 */ 649 #ifdef PETSC_HAVE_SUNDIALS 650 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 651 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 652 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 653 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 654 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 655 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 656 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 657 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 658 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 659 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 660 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 661 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 662 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 663 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 664 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 665 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 666 #endif 667 668 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 669 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 670 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 671 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 672 673 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 674 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 675 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 676 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 677 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 678 679 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 680 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 681 682 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 683 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*); 684 685 #endif 686