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