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