1 /* 2 User interface for the timestepping package. This package 3 is for use in solving time-dependent PDEs. 4 */ 5 #pragma once 6 7 #include <petscsnes.h> 8 #include <petscconvest.h> 9 10 /*I <petscts.h> I*/ 11 12 /* SUBMANSEC = TS */ 13 14 /*S 15 TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 16 17 Level: beginner 18 19 .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()` 20 S*/ 21 typedef struct _p_TS *TS; 22 23 /*J 24 TSType - String with the name of a PETSc `TS` method. 25 26 Level: beginner 27 28 .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()` 29 J*/ 30 typedef const char *TSType; 31 #define TSEULER "euler" 32 #define TSBEULER "beuler" 33 #define TSBASICSYMPLECTIC "basicsymplectic" 34 #define TSPSEUDO "pseudo" 35 #define TSCN "cn" 36 #define TSSUNDIALS "sundials" 37 #define TSRK "rk" 38 #define TSPYTHON "python" 39 #define TSTHETA "theta" 40 #define TSALPHA "alpha" 41 #define TSALPHA2 "alpha2" 42 #define TSGLLE "glle" 43 #define TSGLEE "glee" 44 #define TSSSP "ssp" 45 #define TSARKIMEX "arkimex" 46 #define TSROSW "rosw" 47 #define TSEIMEX "eimex" 48 #define TSMIMEX "mimex" 49 #define TSBDF "bdf" 50 #define TSRADAU5 "radau5" 51 #define TSMPRK "mprk" 52 #define TSDISCGRAD "discgrad" 53 #define TSIRK "irk" 54 #define TSDIRK "dirk" 55 56 /*E 57 TSProblemType - Determines the type of problem this `TS` object is to be used to solve 58 59 Values: 60 + `TS_LINEAR` - a linear ODE or DAE 61 - `TS_NONLINEAR` - a nonlinear ODE or DAE 62 63 Level: beginner 64 65 .seealso: [](ch_ts), `TS`, `TSCreate()` 66 E*/ 67 typedef enum { 68 TS_LINEAR, 69 TS_NONLINEAR 70 } TSProblemType; 71 72 /*E 73 TSEquationType - type of `TS` problem that is solved 74 75 Values: 76 + `TS_EQ_UNSPECIFIED` - (default) 77 . `TS_EQ_EXPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) := M(t) U_t - G(U,t) = 0 78 - `TS_EQ_IMPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0 79 80 Level: beginner 81 82 .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()` 83 E*/ 84 typedef enum { 85 TS_EQ_UNSPECIFIED = -1, 86 TS_EQ_EXPLICIT = 0, 87 TS_EQ_ODE_EXPLICIT = 1, 88 TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 89 TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 90 TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 91 TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 92 TS_EQ_IMPLICIT = 1000, 93 TS_EQ_ODE_IMPLICIT = 1001, 94 TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 95 TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 96 TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 97 TS_EQ_DAE_IMPLICIT_INDEXHI = 1500 98 } TSEquationType; 99 PETSC_EXTERN const char *const *TSEquationTypes; 100 101 /*E 102 TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not 103 104 Values: 105 + `TS_CONVERGED_ITERATING` - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()` 106 . `TS_CONVERGED_TIME` - the final time was reached 107 . `TS_CONVERGED_ITS` - the maximum number of iterations (time-steps) was reached prior to the final time 108 . `TS_CONVERGED_USER` - user requested termination 109 . `TS_CONVERGED_EVENT` - user requested termination on event detection 110 . `TS_CONVERGED_PSEUDO_FATOL` - stops when function norm decreased by a set amount, used only for `TSPSEUDO` 111 . `TS_CONVERGED_PSEUDO_FRTOL` - stops when function norm decreases below a set amount, used only for `TSPSEUDO` 112 . `TS_DIVERGED_NONLINEAR_SOLVE` - too many nonlinear solve failures have occurred 113 . `TS_DIVERGED_STEP_REJECTED` - too many steps were rejected 114 . `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed 115 - `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed 116 117 Level: beginner 118 119 .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()` 120 E*/ 121 typedef enum { 122 TS_CONVERGED_ITERATING = 0, 123 TS_CONVERGED_TIME = 1, 124 TS_CONVERGED_ITS = 2, 125 TS_CONVERGED_USER = 3, 126 TS_CONVERGED_EVENT = 4, 127 TS_CONVERGED_PSEUDO_FATOL = 5, 128 TS_CONVERGED_PSEUDO_FRTOL = 6, 129 TS_DIVERGED_NONLINEAR_SOLVE = -1, 130 TS_DIVERGED_STEP_REJECTED = -2, 131 TSFORWARD_DIVERGED_LINEAR_SOLVE = -3, 132 TSADJOINT_DIVERGED_LINEAR_SOLVE = -4 133 } TSConvergedReason; 134 PETSC_EXTERN const char *const *TSConvergedReasons; 135 136 /*MC 137 TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()` 138 139 Level: beginner 140 141 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()` 142 M*/ 143 144 /*MC 145 TS_CONVERGED_TIME - the final time was reached 146 147 Level: beginner 148 149 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()` 150 M*/ 151 152 /*MC 153 TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time 154 155 Level: beginner 156 157 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()` 158 M*/ 159 160 /*MC 161 TS_CONVERGED_USER - user requested termination 162 163 Level: beginner 164 165 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()` 166 M*/ 167 168 /*MC 169 TS_CONVERGED_EVENT - user requested termination on event detection 170 171 Level: beginner 172 173 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()` 174 M*/ 175 176 /*MC 177 TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO` 178 179 Options Database Key: 180 . -ts_pseudo_frtol <rtol> - use specified rtol 181 182 Level: beginner 183 184 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL` 185 M*/ 186 187 /*MC 188 TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO` 189 190 Options Database Key: 191 . -ts_pseudo_fatol <atol> - use specified atol 192 193 Level: beginner 194 195 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL` 196 M*/ 197 198 /*MC 199 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 200 201 Level: beginner 202 203 Note: 204 See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures. 205 206 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()` 207 M*/ 208 209 /*MC 210 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 211 212 Level: beginner 213 214 Notes: 215 See `TSSetMaxStepRejections()` for how to allow more step rejections. 216 217 .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()` 218 M*/ 219 220 /*E 221 TSExactFinalTimeOption - option for handling of final time step 222 223 Values: 224 + `TS_EXACTFINALTIME_STEPOVER` - Don't do anything if requested final time is exceeded 225 . `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time 226 - `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time requested 227 228 Level: beginner 229 230 .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()` 231 E*/ 232 typedef enum { 233 TS_EXACTFINALTIME_UNSPECIFIED = 0, 234 TS_EXACTFINALTIME_STEPOVER = 1, 235 TS_EXACTFINALTIME_INTERPOLATE = 2, 236 TS_EXACTFINALTIME_MATCHSTEP = 3 237 } TSExactFinalTimeOption; 238 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 239 240 /* Logging support */ 241 PETSC_EXTERN PetscClassId TS_CLASSID; 242 PETSC_EXTERN PetscClassId DMTS_CLASSID; 243 PETSC_EXTERN PetscClassId TSADAPT_CLASSID; 244 245 PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 246 PETSC_EXTERN PetscErrorCode TSFinalizePackage(void); 247 248 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *); 249 PETSC_EXTERN PetscErrorCode TSClone(TS, TS *); 250 PETSC_EXTERN PetscErrorCode TSDestroy(TS *); 251 252 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType); 253 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *); 254 PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec); 255 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **)); 256 PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *)); 257 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 258 259 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]); 260 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]); 261 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]); 262 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 263 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 264 PETSC_EXTERN PetscErrorCode TSReset(TS); 265 266 PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec); 267 PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *); 268 269 PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec); 270 PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *); 271 272 PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *); 273 PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *); 274 PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *); 275 PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec); 276 277 PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *); 278 PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **); 279 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat); 280 PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *); 281 PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool); 282 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *); 283 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *); 284 PETSC_EXTERN PetscErrorCode TSSetIHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *); 285 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 286 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 287 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 288 PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 289 PETSC_EXTERN PetscErrorCode TSSetRHSHessianProduct(TS, Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *, PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *); 290 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 291 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 292 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 293 PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 294 PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec); 295 PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *); 296 PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat); 297 298 /*S 299 TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step) 300 301 Level: advanced 302 303 .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()` 304 S*/ 305 typedef struct _p_TSTrajectory *TSTrajectory; 306 307 /*J 308 TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method 309 310 Level: intermediate 311 312 .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()` 313 J*/ 314 typedef const char *TSTrajectoryType; 315 #define TSTRAJECTORYBASIC "basic" 316 #define TSTRAJECTORYSINGLEFILE "singlefile" 317 #define TSTRAJECTORYMEMORY "memory" 318 #define TSTRAJECTORYVISUALIZATION "visualization" 319 320 PETSC_EXTERN PetscFunctionList TSTrajectoryList; 321 PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID; 322 PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled; 323 324 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS); 325 PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS); 326 PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS); 327 328 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *); 329 PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory); 330 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *); 331 PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer); 332 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType); 333 PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *); 334 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec); 335 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *); 336 PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec); 337 PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *); 338 PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *); 339 PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *); 340 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS); 341 PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS)); 342 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void); 343 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS); 344 PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool); 345 PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool); 346 PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *); 347 PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 348 PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool); 349 PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *); 350 PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool); 351 PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]); 352 PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]); 353 PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *); 354 355 typedef enum { 356 TJ_REVOLVE, 357 TJ_CAMS, 358 TJ_PETSC 359 } TSTrajectoryMemoryType; 360 PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[]; 361 362 PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType); 363 PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt); 364 PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt); 365 PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt); 366 PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt); 367 368 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *); 369 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **); 370 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS() and TSForwardSetSensitivities()", ) PetscErrorCode TSSetCostIntegrand(TS, PetscInt, Vec, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec *, void *), PetscBool, void *); 371 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *); 372 PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec); 373 PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *); 374 PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *); 375 376 PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *); 377 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *); 378 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **)); 379 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS); 380 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *)); 381 382 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *); 383 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat); 384 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *); 385 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *); 386 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS); 387 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt); 388 389 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 390 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 391 PETSC_EXTERN PetscErrorCode TSAdjointReset(TS); 392 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS); 393 PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat); 394 PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS); 395 396 PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat); 397 PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *); 398 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *); 399 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **); 400 PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS); 401 PETSC_EXTERN PetscErrorCode TSForwardReset(TS); 402 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS); 403 PETSC_EXTERN PetscErrorCode TSForwardStep(TS); 404 PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat); 405 PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]); 406 407 PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt); 408 PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *); 409 PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal); 410 PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *); 411 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption); 412 PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *); 413 PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *); 414 PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **); 415 PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **); 416 417 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal); 418 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal); 419 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *); 420 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *); 421 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *); 422 423 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 424 PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 425 426 typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx; 427 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *); 428 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *); 429 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *); 430 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *); 431 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *); 432 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *); 433 434 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *); 435 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *); 436 437 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 438 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *); 439 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *); 440 441 PETSC_EXTERN PetscErrorCode TSStep(TS); 442 PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *); 443 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *); 444 PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec); 445 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *); 446 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType); 447 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *); 448 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason); 449 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *); 450 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *); 451 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *); 452 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *); 453 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt); 454 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *); 455 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt); 456 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool); 457 PETSC_EXTERN PetscErrorCode TSRestartStep(TS); 458 PETSC_EXTERN PetscErrorCode TSRollBack(TS); 459 460 PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]); 461 462 PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *); 463 PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal); 464 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *); 465 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *); 466 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal); 467 PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *); 468 PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt); 469 470 /*S 471 TSRHSFunction - An `TS` right-hand-side evaluation function 472 473 Calling Sequence: 474 + ts - timestep context 475 . t - current time 476 . u - input vector 477 . F - function vector 478 - ctx - [optional] user-defined function context 479 480 Level: beginner 481 482 .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunction()`, 483 `TSIJacobian()`, `TSRHSJacobian()` 484 S*/ 485 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS ts, PetscReal t, Vec u, Vec F, void *ctx); 486 487 /*S 488 TSRHSJacobian - A `TS` right-hand-side Jacobian evaluation function 489 490 Calling Sequence: 491 + ts - the `TS` context obtained from `TSCreate()` 492 . t - current time 493 . u - input vector 494 . Amat - (approximate) Jacobian matrix 495 . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat) 496 - ctx - [optional] user-defined context for matrix evaluation routine 497 498 Level: beginner 499 500 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunction()`, 501 `TSIFunction()`, `TSIJacobian()` 502 S*/ 503 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx); 504 505 /*S 506 TSRHSJacobianP - A function that computes the Jacobian of G w.r.t. the parameters P where U_t 507 = G(U,P,t), as well as the location to store the matrix. 508 509 Calling Sequence: 510 + ts - the `TS` context 511 . t - current timestep 512 . U - input vector (current ODE solution) 513 . A - output matrix 514 - ctx - [optional] user-defined function context 515 516 Level: beginner 517 518 .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()` 519 S*/ 520 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS ts, PetscReal t, Vec U, Mat A, void *ctx); 521 522 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunction, void *); 523 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunction *, void **); 524 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobian, void *); 525 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobian *, void **); 526 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool); 527 528 /*S 529 TSSolutionFunction - A `TS` solution evaluation function 530 531 Calling Sequence: 532 + ts - timestep context 533 . t - current time 534 . u - output vector 535 - ctx - [optional] user-defined function context 536 537 Level: advanced 538 539 .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()` 540 S*/ 541 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS ts, PetscReal t, Vec u, void *ctx); 542 543 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFunction, void *); 544 545 /*S 546 TSForcingFunction - A `TS` forcing function evaluation function 547 548 Calling Sequence: 549 + ts - timestep context 550 . t - current time 551 . f - output vector 552 - ctx - [optional] user-defined function context 553 554 Level: advanced 555 556 .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()` 557 S*/ 558 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS ts, PetscReal t, Vec f, void *ctx); 559 560 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFunction, void *); 561 562 /*S 563 TSIFunction - A `TS` implicit function evaluation function 564 565 Calling Sequence: 566 + ts - the `TS` context obtained from `TSCreate()` 567 . t - time at step/stage being solved 568 . U - state vector 569 . U_t - time derivative of state vector 570 . F - function vector 571 - ctx - [optional] user-defined context for function 572 573 Level: beginner 574 575 .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()` 576 S*/ 577 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS ts, PetscReal t, Vec U, Vec U_tt, Vec F, void *ctx); 578 579 /*S 580 TSIJacobian - A `TS` Jacobian evaluation function 581 582 Calling Sequence: 583 + ts - the `TS` context obtained from `TSCreate()` 584 . t - time at step/stage being solved 585 . U - state vector 586 . U_t - time derivative of state vector 587 . a - shift 588 . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 589 . Pmat - matrix used for constructing preconditioner, usually the same as `Amat` 590 - ctx - [optional] user-defined context for Jacobian evaluation routine 591 592 Level: beginner 593 594 .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunction()`, `TSRHSFunction()`, `TSRHSJacobian()` 595 S*/ 596 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx); 597 598 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunction, void *); 599 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunction *, void **); 600 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobian, void *); 601 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobian *, void **); 602 603 /*S 604 TSI2Function - A `TS` implicit function evaluation function for 2nd order systems 605 606 Calling Sequence: 607 + ts - the `TS` context obtained from `TSCreate()` 608 . t - time at step/stage being solved 609 . U - state vector 610 . U_t - time derivative of state vector 611 . U_tt - second time derivative of state vector 612 . F - function vector 613 - ctx - [optional] user-defined context for matrix evaluation routine (may be `NULL`) 614 615 Level: advanced 616 617 .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()` 618 S*/ 619 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx); 620 621 /*S 622 TSI2Jacobian - A `TS` implicit Jacobian evaluation function for 2nd order systems 623 624 Calling Sequence: 625 + ts - the `TS` context obtained from `TSCreate()` 626 . t - time at step/stage being solved 627 . U - state vector 628 . U_t - time derivative of state vector 629 . U_tt - second time derivative of state vector 630 . v - shift for U_t 631 . a - shift for U_tt 632 . J - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt 633 . jac - matrix from which to construct the preconditioner, may be same as `J` 634 - ctx - [optional] user-defined context for matrix evaluation routine 635 636 Level: advanced 637 638 .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()` 639 S*/ 640 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat Jac, void *ctx); 641 642 PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2Function, void *); 643 PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2Function *, void **); 644 PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2Jacobian, void *); 645 PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2Jacobian *, void **); 646 647 PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS); 648 PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *); 649 PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunction, void *); 650 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *); 651 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]); 652 PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool); 653 PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *); 654 655 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS, PetscReal, Vec, Vec, void *); 656 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS, PetscReal, Vec, Mat, Mat, void *); 657 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *); 658 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 659 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec); 660 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec); 661 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 662 PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat); 663 664 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 665 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal)); 666 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *)); 667 PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS)); 668 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 669 PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *); 670 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 671 PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal); 672 PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *); 673 PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS); 674 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 675 PETSC_EXTERN PetscErrorCode TSResize(TS); 676 PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *); 677 PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec); 678 679 PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec); 680 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec); 681 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *); 682 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *); 683 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *); 684 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal); 685 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *); 686 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *)); 687 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *); 688 689 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *); 690 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *); 691 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *); 692 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal); 693 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *); 694 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *); 695 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *); 696 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal); 697 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 698 699 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]); 700 PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]); 701 702 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec); 703 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat); 704 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool); 705 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool); 706 PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec); 707 PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat); 708 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *); 709 710 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec); 711 712 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 713 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunction, void *); 714 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunction *, void **); 715 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *)); 716 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobian, void *); 717 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobian *, void **); 718 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *)); 719 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunction, void *); 720 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunction *, void **); 721 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *)); 722 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobian, void *); 723 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobian *, void **); 724 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *)); 725 PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2Function, void *); 726 PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2Function *, void **); 727 PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *)); 728 PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2Jacobian, void *); 729 PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2Jacobian *, void **); 730 PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *)); 731 732 /*S 733 TSTransientVariable - A function to transform from state to transient variables 734 735 Calling Sequence: 736 + ts - timestep context 737 . p - input vector (primitive form) 738 . c - output vector, transient variables (conservative form) 739 - ctx - [optional] user-defined function context 740 741 Level: advanced 742 743 .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()` 744 S*/ 745 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS ts, Vec p, Vec c, void *ctx); 746 747 PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariable, void *); 748 PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariable, void *); 749 PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariable *, void *); 750 PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec); 751 PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *); 752 753 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFunction, void *); 754 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFunction *, void **); 755 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFunction, void *); 756 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFunction *, void **); 757 PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *); 758 PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *); 759 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec); 760 761 PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **); 762 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *); 763 PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **); 764 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *); 765 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **); 766 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 767 PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM); 768 PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM); 769 PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM); 770 771 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer)); 772 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer)); 773 774 /*S 775 DMDATSRHSFunctionLocal - A local `TS` right hand side residual evaluation function for use with `DMDA` 776 777 Calling Sequence: 778 + info - defines the subdomain to evaluate the residual on 779 . t - time at which to evaluate residual 780 . x - array of local state information 781 . f - output array of local residual information 782 - ctx - optional user context 783 784 Level: beginner 785 786 .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunction()`, `DMDATSRHSJacobianLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()` 787 S*/ 788 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx); 789 790 /*S 791 DMDATSRHSJacobianLocal - A local residual evaluation function for use with `DMDA` 792 793 Calling Sequence: 794 + info - defines the subdomain to evaluate the residual on 795 . t - time at which to evaluate residual 796 . x - array of local state information 797 . J - Jacobian matrix 798 . B - matrix from which to construct the preconditioner; often same as `J` 799 - ctx - optional context 800 801 Level: beginner 802 803 .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobian()`, `DMDATSRHSFunctionLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()` 804 S*/ 805 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx); 806 807 /*S 808 DMDATSIFunctionLocal - A local residual evaluation function for use with `DMDA` 809 810 Calling Sequence: 811 + info - defines the subdomain to evaluate the residual on 812 . t - time at which to evaluate residual 813 . x - array of local state information 814 . xdot - array of local time derivative information 815 . imode - output array of local function evaluation information 816 - ctx - optional context 817 818 Level: beginner 819 820 .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocal()`, `TSIFunction()` 821 S*/ 822 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx); 823 824 /*S 825 DMDATSIJacobianLocal - A local residual evaluation function for use with `DMDA` 826 827 Calling Sequence: 828 + info - defines the subdomain to evaluate the residual on 829 . t - time at which to evaluate the jacobian 830 . x - array of local state information 831 . xdot - time derivative at this state 832 . shift - see `TSSetIJacobian()` for the meaning of this parameter 833 . J - Jacobian matrix 834 . B - matrix from which to construct the preconditioner; often same as `J` 835 - ctx - optional context 836 837 Level: beginner 838 839 .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobian()`, `DMDATSIFunctionLocal()`, `DMDATSRHSFunctionLocal()`, `DMDATSRHSJacob ianlocal()` 840 S*/ 841 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx); 842 843 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocal, void *); 844 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocal, void *); 845 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocal, void *); 846 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocal, void *); 847 848 typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx; 849 typedef struct { 850 Vec ray; 851 VecScatter scatter; 852 PetscViewer viewer; 853 TSMonitorLGCtx lgctx; 854 } TSMonitorDMDARayCtx; 855 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **); 856 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *); 857 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *); 858 859 /* Dynamic creation and loading functions */ 860 PETSC_EXTERN PetscFunctionList TSList; 861 PETSC_EXTERN PetscErrorCode TSGetType(TS, TSType *); 862 PETSC_EXTERN PetscErrorCode TSSetType(TS, TSType); 863 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 864 865 PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *); 866 PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES); 867 PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *); 868 869 PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer); 870 PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer); 871 PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]); 872 PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]); 873 874 #define TS_FILE_CLASSID 1211225 875 876 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *); 877 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *); 878 879 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *); 880 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *); 881 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *); 882 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *); 883 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *); 884 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **); 885 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *); 886 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *); 887 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *); 888 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 889 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 890 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *); 891 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *); 892 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *); 893 PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 894 PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 895 896 struct _n_TSMonitorLGCtxNetwork { 897 PetscInt nlg; 898 PetscDrawLG *lg; 899 PetscBool semilogy; 900 PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */ 901 }; 902 typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork; 903 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *); 904 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *); 905 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *); 906 907 typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx; 908 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *); 909 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *); 910 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *); 911 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *); 912 913 typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx; 914 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *); 915 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *); 916 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *); 917 918 typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx; 919 PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *); 920 PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *); 921 PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 922 923 typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx; 924 PETSC_EXTERN PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *); 925 PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 926 PETSC_EXTERN PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *); 927 PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 928 929 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *); 930 PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal); 931 PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal); 932 PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt) 933 { 934 return TSSetPostEventSecondStep(ts, dt); 935 } 936 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]); 937 PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *); 938 939 /*J 940 TSSSPType - string with the name of a `TSSSP` scheme. 941 942 Level: beginner 943 944 .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP` 945 J*/ 946 typedef const char *TSSSPType; 947 #define TSSSPRKS2 "rks2" 948 #define TSSSPRKS3 "rks3" 949 #define TSSSPRK104 "rk104" 950 951 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS, TSSSPType); 952 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS, TSSSPType *); 953 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS, PetscInt); 954 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS, PetscInt *); 955 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 956 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 957 PETSC_EXTERN PetscFunctionList TSSSPList; 958 959 /*S 960 TSAdapt - Abstract object that manages time-step adaptivity 961 962 Level: beginner 963 964 .seealso: [](ch_ts), `TS`, `TSAdaptCreate()`, `TSAdaptType` 965 S*/ 966 typedef struct _p_TSAdapt *TSAdapt; 967 968 /*J 969 TSAdaptType - String with the name of `TSAdapt` scheme. 970 971 Level: beginner 972 973 .seealso: [](ch_ts), `TSAdaptSetType()`, `TS`, `TSAdapt` 974 J*/ 975 typedef const char *TSAdaptType; 976 #define TSADAPTNONE "none" 977 #define TSADAPTBASIC "basic" 978 #define TSADAPTDSP "dsp" 979 #define TSADAPTCFL "cfl" 980 #define TSADAPTGLEE "glee" 981 #define TSADAPTHISTORY "history" 982 983 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *); 984 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt)); 985 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 986 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 987 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *); 988 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType); 989 PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *); 990 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]); 991 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 992 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool); 993 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **); 994 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *); 995 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *); 996 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer); 997 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer); 998 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *); 999 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt); 1000 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *); 1001 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool); 1002 PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool); 1003 PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal); 1004 PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *); 1005 PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal); 1006 PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *); 1007 PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal); 1008 PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *); 1009 PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal); 1010 PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *); 1011 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal); 1012 PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *); 1013 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *)); 1014 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool); 1015 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool); 1016 PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *); 1017 PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt); 1018 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *); 1019 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal); 1020 1021 /*S 1022 TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE` 1023 1024 Level: beginner 1025 1026 Developer Note: 1027 This functionality should be replaced by the `TSAdapt`. 1028 1029 .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType` 1030 S*/ 1031 typedef struct _p_TSGLLEAdapt *TSGLLEAdapt; 1032 1033 /*J 1034 TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme 1035 1036 Level: beginner 1037 1038 Developer Note: 1039 This functionality should be replaced by the `TSAdaptType`. 1040 1041 .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS` 1042 J*/ 1043 typedef const char *TSGLLEAdaptType; 1044 #define TSGLLEADAPT_NONE "none" 1045 #define TSGLLEADAPT_SIZE "size" 1046 #define TSGLLEADAPT_BOTH "both" 1047 1048 PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt)); 1049 PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void); 1050 PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void); 1051 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *); 1052 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType); 1053 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]); 1054 PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *); 1055 PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer); 1056 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *); 1057 PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *); 1058 1059 /*J 1060 TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme 1061 1062 Level: beginner 1063 1064 .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept` 1065 J*/ 1066 typedef const char *TSGLLEAcceptType; 1067 #define TSGLLEACCEPT_ALWAYS "always" 1068 1069 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS, PetscReal, PetscReal, const PetscReal[], PetscBool *); 1070 PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFunction); 1071 1072 /*J 1073 TSGLLEType - string with the name of a General Linear `TSGLLE` type 1074 1075 Level: beginner 1076 1077 .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept` 1078 J*/ 1079 typedef const char *TSGLLEType; 1080 #define TSGLLE_IRKS "irks" 1081 1082 PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS)); 1083 PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void); 1084 PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void); 1085 PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType); 1086 PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *); 1087 PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType); 1088 1089 /*J 1090 TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type 1091 1092 Level: beginner 1093 1094 .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()` 1095 J*/ 1096 #define TSEIMEXType char * 1097 1098 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt); 1099 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt); 1100 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool); 1101 1102 /*J 1103 TSRKType - String with the name of a Runge-Kutta `TSRK` type 1104 1105 Level: beginner 1106 1107 .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()` 1108 J*/ 1109 typedef const char *TSRKType; 1110 #define TSRK1FE "1fe" 1111 #define TSRK2A "2a" 1112 #define TSRK2B "2b" 1113 #define TSRK3 "3" 1114 #define TSRK3BS "3bs" 1115 #define TSRK4 "4" 1116 #define TSRK5F "5f" 1117 #define TSRK5DP "5dp" 1118 #define TSRK5BS "5bs" 1119 #define TSRK6VR "6vr" 1120 #define TSRK7VR "7vr" 1121 #define TSRK8VR "8vr" 1122 1123 PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *); 1124 PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *); 1125 PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType); 1126 PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *); 1127 PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool); 1128 PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *); 1129 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1130 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 1131 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 1132 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 1133 1134 /*J 1135 TSMPRKType - String with the name of a Partitioned Runge-Kutta `TSMPRK` type 1136 1137 Level: beginner 1138 1139 .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()` 1140 J*/ 1141 typedef const char *TSMPRKType; 1142 #define TSMPRK2A22 "2a22" 1143 #define TSMPRK2A23 "2a23" 1144 #define TSMPRK2A32 "2a32" 1145 #define TSMPRK2A33 "2a33" 1146 #define TSMPRKP2 "p2" 1147 #define TSMPRKP3 "p3" 1148 1149 PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *); 1150 PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType); 1151 PETSC_EXTERN PetscErrorCode TSMPRKRegister(TSMPRKType, PetscInt, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscInt[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscInt[], const PetscReal[], const PetscReal[], const PetscReal[]); 1152 PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void); 1153 PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void); 1154 PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void); 1155 1156 /*J 1157 TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type 1158 1159 Level: beginner 1160 1161 .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()` 1162 J*/ 1163 typedef const char *TSIRKType; 1164 #define TSIRKGAUSS "gauss" 1165 1166 PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *); 1167 PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType); 1168 PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *); 1169 PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt); 1170 PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS)); 1171 PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *); 1172 PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void); 1173 PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void); 1174 PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void); 1175 1176 /*J 1177 TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type 1178 1179 Level: beginner 1180 1181 .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()` 1182 J*/ 1183 typedef const char *TSGLEEType; 1184 #define TSGLEEi1 "BE1" 1185 #define TSGLEE23 "23" 1186 #define TSGLEE24 "24" 1187 #define TSGLEE25I "25i" 1188 #define TSGLEE35 "35" 1189 #define TSGLEEEXRK2A "exrk2a" 1190 #define TSGLEERK32G1 "rk32g1" 1191 #define TSGLEERK285EX "rk285ex" 1192 1193 /*J 1194 TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type 1195 1196 Level: beginner 1197 1198 .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()` 1199 J*/ 1200 PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *); 1201 PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType); 1202 PETSC_EXTERN PetscErrorCode TSGLEERegister(TSGLEEType, PetscInt, PetscInt, PetscInt, PetscReal, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1203 PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void); 1204 PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void); 1205 PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void); 1206 PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void); 1207 1208 /*J 1209 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type 1210 1211 Level: beginner 1212 1213 .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()` 1214 J*/ 1215 typedef const char *TSARKIMEXType; 1216 #define TSARKIMEX1BEE "1bee" 1217 #define TSARKIMEXA2 "a2" 1218 #define TSARKIMEXL2 "l2" 1219 #define TSARKIMEXARS122 "ars122" 1220 #define TSARKIMEX2C "2c" 1221 #define TSARKIMEX2D "2d" 1222 #define TSARKIMEX2E "2e" 1223 #define TSARKIMEXPRSSP2 "prssp2" 1224 #define TSARKIMEX3 "3" 1225 #define TSARKIMEXBPR3 "bpr3" 1226 #define TSARKIMEXARS443 "ars443" 1227 #define TSARKIMEX4 "4" 1228 #define TSARKIMEX5 "5" 1229 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *); 1230 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType); 1231 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool); 1232 PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *); 1233 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[]); 1234 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 1235 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 1236 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 1237 1238 /*J 1239 TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type 1240 1241 Level: beginner 1242 1243 .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()` 1244 J*/ 1245 typedef const char *TSDIRKType; 1246 #define TSDIRKS212 "s212" 1247 #define TSDIRKES122SAL "es122sal" 1248 #define TSDIRKES213SAL "es213sal" 1249 #define TSDIRKES324SAL "es324sal" 1250 #define TSDIRKES325SAL "es325sal" 1251 #define TSDIRK657A "657a" 1252 #define TSDIRKES648SA "es648sa" 1253 #define TSDIRK658A "658a" 1254 #define TSDIRKS659A "s659a" 1255 #define TSDIRK7510SAL "7510sal" 1256 #define TSDIRKES7510SA "es7510sa" 1257 #define TSDIRK759A "759a" 1258 #define TSDIRKS7511SAL "s7511sal" 1259 #define TSDIRK8614A "8614a" 1260 #define TSDIRK8616SAL "8616sal" 1261 #define TSDIRKES8516SAL "es8516sal" 1262 PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *); 1263 PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType); 1264 PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1265 1266 /*J 1267 TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type 1268 1269 Level: beginner 1270 1271 .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()` 1272 J*/ 1273 typedef const char *TSRosWType; 1274 #define TSROSW2M "2m" 1275 #define TSROSW2P "2p" 1276 #define TSROSWRA3PW "ra3pw" 1277 #define TSROSWRA34PW2 "ra34pw2" 1278 #define TSROSWRODAS3 "rodas3" 1279 #define TSROSWSANDU3 "sandu3" 1280 #define TSROSWASSP3P3S1C "assp3p3s1c" 1281 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 1282 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 1283 #define TSROSWARK3 "ark3" 1284 #define TSROSWTHETA1 "theta1" 1285 #define TSROSWTHETA2 "theta2" 1286 #define TSROSWGRK4T "grk4t" 1287 #define TSROSWSHAMP4 "shamp4" 1288 #define TSROSWVELDD4 "veldd4" 1289 #define TSROSW4L "4l" 1290 1291 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *); 1292 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType); 1293 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool); 1294 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1295 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal); 1296 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 1297 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 1298 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 1299 1300 PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt); 1301 PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *); 1302 1303 /*J 1304 TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type 1305 1306 Level: beginner 1307 1308 .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()` 1309 J*/ 1310 typedef const char *TSBasicSymplecticType; 1311 #define TSBASICSYMPLECTICSIEULER "1" 1312 #define TSBASICSYMPLECTICVELVERLET "2" 1313 #define TSBASICSYMPLECTIC3 "3" 1314 #define TSBASICSYMPLECTIC4 "4" 1315 PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType); 1316 PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *); 1317 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]); 1318 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void); 1319 PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void); 1320 PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void); 1321 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void); 1322 1323 /*J 1324 TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy), 1325 but also has the property for some systems of monotonicity in a functional. 1326 1327 Level: beginner 1328 1329 .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()` 1330 J*/ 1331 PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *); 1332 PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *); 1333 PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *); 1334 PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool); 1335 1336 /* 1337 PETSc interface to Sundials 1338 */ 1339 #ifdef PETSC_HAVE_SUNDIALS2 1340 typedef enum { 1341 SUNDIALS_ADAMS = 1, 1342 SUNDIALS_BDF = 2 1343 } TSSundialsLmmType; 1344 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 1345 typedef enum { 1346 SUNDIALS_MODIFIED_GS = 1, 1347 SUNDIALS_CLASSICAL_GS = 2 1348 } TSSundialsGramSchmidtType; 1349 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 1350 1351 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType); 1352 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *); 1353 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal); 1354 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal); 1355 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal); 1356 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *); 1357 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType); 1358 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt); 1359 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal); 1360 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool); 1361 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt); 1362 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt); 1363 PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool); 1364 #endif 1365 1366 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal); 1367 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *); 1368 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *); 1369 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool); 1370 1371 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal); 1372 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal); 1373 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *); 1374 1375 PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal); 1376 PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal); 1377 PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 1378 1379 /*S 1380 TSAlpha2Predictor - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in 1381 a second-order generalized-alpha time integrator. 1382 1383 Calling Sequence: 1384 + ts - the `TS` context obtained from `TSCreate()` 1385 . X0 - the previous time step's state vector 1386 . V0 - the previous time step's first derivative of the state vector 1387 . A0 - the previous time step's second derivative of the state vector 1388 . X1 - the vector into which the initial guess for the current time step will be written 1389 - ctx - [optional] user-defined context for the predictor evaluation routine (may be `NULL`) 1390 1391 Level: intermediate 1392 1393 .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()` 1394 S*/ 1395 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSAlpha2Predictor)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx); 1396 PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2Predictor, void *ctx); 1397 1398 PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM); 1399 PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *); 1400 1401 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *); 1402 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *); 1403 1404 PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *); 1405 PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *); 1406 1407 PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec)); 1408 PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec)); 1409 PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec); 1410 PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec)); 1411 PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec)); 1412 PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec); 1413 PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool); 1414 1415 PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure); 1416