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