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