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