1818ad0c1SBarry Smith /* 2316643e7SJed Brown User interface for the timestepping package. This package 3f64a0f93SLois Curfman McInnes is for use in solving time-dependent PDEs. 4818ad0c1SBarry Smith */ 5a4963045SJacob Faibussowitsch #pragma once 6ac09b921SBarry Smith 72c8e378dSBarry Smith #include <petscsnes.h> 8900f6b5bSMatthew G. Knepley #include <petscconvest.h> 9818ad0c1SBarry Smith 1014d0ab18SJacob Faibussowitsch /*I <petscts.h> I*/ 1114d0ab18SJacob Faibussowitsch 12ac09b921SBarry Smith /* SUBMANSEC = TS */ 13ac09b921SBarry Smith 14435da068SBarry Smith /*S 150b4b7b1cSBarry Smith TS - Abstract PETSc object that manages integrating an ODE. 16435da068SBarry Smith 17435da068SBarry Smith Level: beginner 18435da068SBarry Smith 191cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSCreate()`, `TSSetType()`, `TSType`, `SNES`, `KSP`, `PC`, `TSDestroy()` 20435da068SBarry Smith S*/ 21f09e8eb9SSatish Balay typedef struct _p_TS *TS; 22435da068SBarry Smith 2376bdecfbSBarry Smith /*J 240b4b7b1cSBarry Smith TSType - String with the name of a PETSc `TS` method. These are all the time/ODE integrators that PETSc provides. 25435da068SBarry Smith 26435da068SBarry Smith Level: beginner 27435da068SBarry Smith 280b4b7b1cSBarry Smith Note: 290b4b7b1cSBarry Smith Use `TSSetType()` or the options database key `-ts_type` to set the ODE integrator method to use with a given `TS` object 300b4b7b1cSBarry Smith 311cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()` 3276bdecfbSBarry Smith J*/ 3319fd82e9SBarry Smith typedef const char *TSType; 349596e0b4SJed Brown #define TSEULER "euler" 359596e0b4SJed Brown #define TSBEULER "beuler" 36e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic" 379596e0b4SJed Brown #define TSPSEUDO "pseudo" 384d91e141SJed Brown #define TSCN "cn" 399596e0b4SJed Brown #define TSSUNDIALS "sundials" 404d91e141SJed Brown #define TSRK "rk" 419596e0b4SJed Brown #define TSPYTHON "python" 429596e0b4SJed Brown #define TSTHETA "theta" 4388df8a41SLisandro Dalcin #define TSALPHA "alpha" 44818efac9SLisandro Dalcin #define TSALPHA2 "alpha2" 4526d28e4eSEmil Constantinescu #define TSGLLE "glle" 46b6a60446SDebojyoti Ghosh #define TSGLEE "glee" 47ef7bb5aaSJed Brown #define TSSSP "ssp" 488a381b04SJed Brown #define TSARKIMEX "arkimex" 49e27a552bSJed Brown #define TSROSW "rosw" 507d5697caSHong Zhang #define TSEIMEX "eimex" 51abadfa0fSMatthew G. Knepley #define TSMIMEX "mimex" 52211a84d6SLisandro Dalcin #define TSBDF "bdf" 53d249a78fSBarry Smith #define TSRADAU5 "radau5" 544b84eec9SHong Zhang #define TSMPRK "mprk" 5540be0ff1SMatthew G. Knepley #define TSDISCGRAD "discgrad" 56d2567f34SHong Zhang #define TSIRK "irk" 57d27334e2SStefano Zampini #define TSDIRK "dirk" 5840be0ff1SMatthew G. Knepley 59435da068SBarry Smith /*E 6087497f52SBarry Smith TSProblemType - Determines the type of problem this `TS` object is to be used to solve 61435da068SBarry Smith 62195e9b02SBarry Smith Values: 63195e9b02SBarry Smith + `TS_LINEAR` - a linear ODE or DAE 64195e9b02SBarry Smith - `TS_NONLINEAR` - a nonlinear ODE or DAE 65195e9b02SBarry Smith 66435da068SBarry Smith Level: beginner 67435da068SBarry Smith 681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()` 69435da068SBarry Smith E*/ 709371c9d4SSatish Balay typedef enum { 719371c9d4SSatish Balay TS_LINEAR, 729371c9d4SSatish Balay TS_NONLINEAR 739371c9d4SSatish Balay } TSProblemType; 74818ad0c1SBarry Smith 75b93fea0eSJed Brown /*E 7687497f52SBarry Smith TSEquationType - type of `TS` problem that is solved 77e817cc15SEmil Constantinescu 78195e9b02SBarry Smith Values: 79195e9b02SBarry Smith + `TS_EQ_UNSPECIFIED` - (default) 8016a05f60SBarry Smith . `TS_EQ_EXPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) := M(t) U_t - G(U,t) = 0 8116a05f60SBarry Smith - `TS_EQ_IMPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0 82e817cc15SEmil Constantinescu 8395bd0b28SBarry Smith Level: beginner 8495bd0b28SBarry Smith 851cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()` 86e817cc15SEmil Constantinescu E*/ 87e817cc15SEmil Constantinescu typedef enum { 88e817cc15SEmil Constantinescu TS_EQ_UNSPECIFIED = -1, 89e817cc15SEmil Constantinescu TS_EQ_EXPLICIT = 0, 90e817cc15SEmil Constantinescu TS_EQ_ODE_EXPLICIT = 1, 91e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 92e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 93e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 94e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 95e817cc15SEmil Constantinescu TS_EQ_IMPLICIT = 1000, 96e817cc15SEmil Constantinescu TS_EQ_ODE_IMPLICIT = 1001, 97e817cc15SEmil Constantinescu TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 98e817cc15SEmil Constantinescu TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 99e817cc15SEmil Constantinescu TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 10019436ca2SJed Brown TS_EQ_DAE_IMPLICIT_INDEXHI = 1500 101e817cc15SEmil Constantinescu } TSEquationType; 102e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes; 103e817cc15SEmil Constantinescu 104e817cc15SEmil Constantinescu /*E 10516a05f60SBarry Smith TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not 106b93fea0eSJed Brown 107195e9b02SBarry Smith Values: 108195e9b02SBarry Smith + `TS_CONVERGED_ITERATING` - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()` 109195e9b02SBarry Smith . `TS_CONVERGED_TIME` - the final time was reached 110195e9b02SBarry Smith . `TS_CONVERGED_ITS` - the maximum number of iterations (time-steps) was reached prior to the final time 111195e9b02SBarry Smith . `TS_CONVERGED_USER` - user requested termination 112195e9b02SBarry Smith . `TS_CONVERGED_EVENT` - user requested termination on event detection 113195e9b02SBarry Smith . `TS_CONVERGED_PSEUDO_FATOL` - stops when function norm decreased by a set amount, used only for `TSPSEUDO` 114195e9b02SBarry Smith . `TS_CONVERGED_PSEUDO_FRTOL` - stops when function norm decreases below a set amount, used only for `TSPSEUDO` 115195e9b02SBarry Smith . `TS_DIVERGED_NONLINEAR_SOLVE` - too many nonlinear solve failures have occurred 116195e9b02SBarry Smith . `TS_DIVERGED_STEP_REJECTED` - too many steps were rejected 117195e9b02SBarry Smith . `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed 118195e9b02SBarry Smith - `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed 119195e9b02SBarry Smith 120b93fea0eSJed Brown Level: beginner 121b93fea0eSJed Brown 1221cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()` 123b93fea0eSJed Brown E*/ 124193ac0bcSJed Brown typedef enum { 125193ac0bcSJed Brown TS_CONVERGED_ITERATING = 0, 126193ac0bcSJed Brown TS_CONVERGED_TIME = 1, 127193ac0bcSJed Brown TS_CONVERGED_ITS = 2, 128d6ad946cSShri Abhyankar TS_CONVERGED_USER = 3, 1292b7db910SShri Abhyankar TS_CONVERGED_EVENT = 4, 1303118ae5eSBarry Smith TS_CONVERGED_PSEUDO_FATOL = 5, 1313118ae5eSBarry Smith TS_CONVERGED_PSEUDO_FRTOL = 6, 132193ac0bcSJed Brown TS_DIVERGED_NONLINEAR_SOLVE = -1, 133b5b4017aSHong Zhang TS_DIVERGED_STEP_REJECTED = -2, 13405c9950eSHong Zhang TSFORWARD_DIVERGED_LINEAR_SOLVE = -3, 13505c9950eSHong Zhang TSADJOINT_DIVERGED_LINEAR_SOLVE = -4 136193ac0bcSJed Brown } TSConvergedReason; 137014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons; 1381957e957SBarry Smith 139b93fea0eSJed Brown /*MC 14087497f52SBarry Smith TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()` 141b93fea0eSJed Brown 142b93fea0eSJed Brown Level: beginner 143b93fea0eSJed Brown 1441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()` 145b93fea0eSJed Brown M*/ 146b93fea0eSJed Brown 147b93fea0eSJed Brown /*MC 148b93fea0eSJed Brown TS_CONVERGED_TIME - the final time was reached 149b93fea0eSJed Brown 150b93fea0eSJed Brown Level: beginner 151b93fea0eSJed Brown 1521cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()` 153b93fea0eSJed Brown M*/ 154b93fea0eSJed Brown 155b93fea0eSJed Brown /*MC 1561957e957SBarry Smith TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time 157b93fea0eSJed Brown 158b93fea0eSJed Brown Level: beginner 159b93fea0eSJed Brown 1601cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()` 161b93fea0eSJed Brown M*/ 1621957e957SBarry Smith 163d6ad946cSShri Abhyankar /*MC 164d6ad946cSShri Abhyankar TS_CONVERGED_USER - user requested termination 165d6ad946cSShri Abhyankar 166d6ad946cSShri Abhyankar Level: beginner 167d6ad946cSShri Abhyankar 1681cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()` 169b93fea0eSJed Brown M*/ 170b93fea0eSJed Brown 171b93fea0eSJed Brown /*MC 172d76fc68cSShri Abhyankar TS_CONVERGED_EVENT - user requested termination on event detection 173d76fc68cSShri Abhyankar 174d76fc68cSShri Abhyankar Level: beginner 175d76fc68cSShri Abhyankar 1761cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()` 177d76fc68cSShri Abhyankar M*/ 178d76fc68cSShri Abhyankar 179d76fc68cSShri Abhyankar /*MC 18087497f52SBarry Smith TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO` 1813118ae5eSBarry Smith 1823c7db156SBarry Smith Options Database Key: 183d5b43468SJose E. Roman . -ts_pseudo_frtol <rtol> - use specified rtol 1843118ae5eSBarry Smith 185195e9b02SBarry Smith Level: beginner 186195e9b02SBarry Smith 1871cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL` 1883118ae5eSBarry Smith M*/ 1893118ae5eSBarry Smith 1903118ae5eSBarry Smith /*MC 19187497f52SBarry Smith TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO` 1923118ae5eSBarry Smith 1933c7db156SBarry Smith Options Database Key: 19467b8a455SSatish Balay . -ts_pseudo_fatol <atol> - use specified atol 1953118ae5eSBarry Smith 196195e9b02SBarry Smith Level: beginner 197195e9b02SBarry Smith 1981cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL` 1993118ae5eSBarry Smith M*/ 2003118ae5eSBarry Smith 2013118ae5eSBarry Smith /*MC 202b93fea0eSJed Brown TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 203b93fea0eSJed Brown 204b93fea0eSJed Brown Level: beginner 205b93fea0eSJed Brown 206195e9b02SBarry Smith Note: 207195e9b02SBarry Smith See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures. 2081957e957SBarry Smith 2091cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()` 210b93fea0eSJed Brown M*/ 211b93fea0eSJed Brown 212b93fea0eSJed Brown /*MC 213b93fea0eSJed Brown TS_DIVERGED_STEP_REJECTED - too many steps were rejected 214b93fea0eSJed Brown 215b93fea0eSJed Brown Level: beginner 216b93fea0eSJed Brown 21795452b02SPatrick Sanan Notes: 218195e9b02SBarry Smith See `TSSetMaxStepRejections()` for how to allow more step rejections. 2191957e957SBarry Smith 2201cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()` 221b93fea0eSJed Brown M*/ 222b93fea0eSJed Brown 223d6ad946cSShri Abhyankar /*E 224d6ad946cSShri Abhyankar TSExactFinalTimeOption - option for handling of final time step 225d6ad946cSShri Abhyankar 226195e9b02SBarry Smith Values: 22716a05f60SBarry Smith + `TS_EXACTFINALTIME_STEPOVER` - Don't do anything if requested final time is exceeded 228195e9b02SBarry Smith . `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time 22916a05f60SBarry Smith - `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time requested 230195e9b02SBarry Smith 231d6ad946cSShri Abhyankar Level: beginner 232d6ad946cSShri Abhyankar 2331cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()` 234d6ad946cSShri Abhyankar E*/ 2359371c9d4SSatish Balay typedef enum { 2369371c9d4SSatish Balay TS_EXACTFINALTIME_UNSPECIFIED = 0, 2379371c9d4SSatish Balay TS_EXACTFINALTIME_STEPOVER = 1, 2389371c9d4SSatish Balay TS_EXACTFINALTIME_INTERPOLATE = 2, 2399371c9d4SSatish Balay TS_EXACTFINALTIME_MATCHSTEP = 3 2409371c9d4SSatish Balay } TSExactFinalTimeOption; 241d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 242d6ad946cSShri Abhyankar 243000e7ae3SMatthew Knepley /* Logging support */ 244014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID; 245d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID; 2461b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID; 247000e7ae3SMatthew Knepley 248607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 249560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void); 2508ba1e511SMatthew Knepley 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *); 252baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *); 254818ad0c1SBarry Smith 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec); 25849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscCtxDestroyFn *); 259721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *)); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 261818ad0c1SBarry Smith 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS); 268818ad0c1SBarry Smith 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *); 271818ad0c1SBarry Smith 272efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec); 273efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *); 27457df6a1bSDebojyoti Ghosh 275642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *); 276642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *); 2770a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *); 27857df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec); 2796e2e232bSDebojyoti Ghosh 280a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *); 281cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **); 282a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat); 28333f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *); 28433f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool); 285edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *); 286edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *); 2879371c9d4SSatish Balay 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 *); 288b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 289b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 290b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 291b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 2929371c9d4SSatish Balay 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 *); 293b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 294b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 295b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 296b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 297b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec); 298b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *); 299ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat); 300a05bf03eSHong Zhang 301bc952696SBarry Smith /*S 302e91eccc2SStefano Zampini TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step) 303bc952696SBarry Smith 304bc952696SBarry Smith Level: advanced 305bc952696SBarry Smith 3061cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()` 307bc952696SBarry Smith S*/ 308bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory; 309bc952696SBarry Smith 310bc952696SBarry Smith /*J 311195e9b02SBarry Smith TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method 312bc952696SBarry Smith 313bc952696SBarry Smith Level: intermediate 314bc952696SBarry Smith 3151cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()` 316bc952696SBarry Smith J*/ 317bc952696SBarry Smith typedef const char *TSTrajectoryType; 318bc952696SBarry Smith #define TSTRAJECTORYBASIC "basic" 3191c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE "singlefile" 3209a53571cSHong Zhang #define TSTRAJECTORYMEMORY "memory" 3212b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization" 322bc952696SBarry Smith 323bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList; 324bc952696SBarry Smith PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID; 325bc952696SBarry Smith PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled; 326bc952696SBarry Smith 327bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS); 3282d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS); 32967a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS); 330bc952696SBarry Smith 331bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *); 3329a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory); 333bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *); 334560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer); 335fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType); 336881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *); 337bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec); 338c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *); 339fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec); 340fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *); 341fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *); 342fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *); 343972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS); 344560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS)); 345bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void); 34668bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS); 3479ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool); 3482bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool); 34978fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *); 35008347785SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 351fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool); 352fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *); 35364fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool); 35464e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]); 35564e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]); 35678fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *); 357bc952696SBarry Smith 3584bf303faSJacob Faibussowitsch typedef enum { 3594bf303faSJacob Faibussowitsch TJ_REVOLVE, 3604bf303faSJacob Faibussowitsch TJ_CAMS, 3614bf303faSJacob Faibussowitsch TJ_PETSC 3624bf303faSJacob Faibussowitsch } TSTrajectoryMemoryType; 3634bf303faSJacob Faibussowitsch PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[]; 3644bf303faSJacob Faibussowitsch 3654bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType); 3664bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt); 3674bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt); 3684bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt); 3694bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt); 3704bf303faSJacob Faibussowitsch 371dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *); 372dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **); 373edd03b47SJacob Faibussowitsch 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 *); 374dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *); 375022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec); 376cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *); 377cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *); 378d4aa0a58SBarry Smith 379dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *); 380a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *); 38149abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscCtxDestroyFn *); 382a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS); 383a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *)); 384a05bf03eSHong Zhang 385edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *); 386edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat); 387edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *); 388edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *); 3892c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS); 39073b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt); 391d4aa0a58SBarry Smith 39273b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 393f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 394ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS); 395999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS); 396ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat); 397ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS); 398715f1b00SHong Zhang 39913b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat); 40013b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *); 401edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *); 402edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **); 403715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS); 4047adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS); 405999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS); 406715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS); 407b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat); 408cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]); 409c235aa8dSHong Zhang 410618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt); 411618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *); 412618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal); 413618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *); 41449354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption); 415f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *); 4164a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *); 4174a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **); 4184a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **); 419*8343f784SJames Wright PETSC_EXTERN PetscErrorCode TSSetEvaluationTimes(TS, PetscInt, PetscReal *); 420*8343f784SJames Wright PETSC_EXTERN PetscErrorCode TSGetEvaluationTimes(TS, PetscInt *, const PetscReal **); 421*8343f784SJames Wright PETSC_EXTERN PetscErrorCode TSGetEvaluationTimesSolutions(TS, PetscInt *, const PetscReal **, Vec **); 422818ad0c1SBarry Smith 423edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal); 424edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal); 425edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *); 426edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *); 427edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *); 42819eac22cSLisandro Dalcin 429721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 430cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 43183a4ac43SBarry Smith 43283a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx; 43383a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *); 43483a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *); 43583a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *); 4362d5ee99bSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *); 43783a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *); 4380ed3bfb6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *); 43983a4ac43SBarry Smith 440721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *); 4419110b2e7SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *); 4429110b2e7SHong Zhang 443721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 4447f27e910SStefano Zampini 4457f27e910SStefano Zampini typedef struct _n_TSMonitorVTKCtx *TSMonitorVTKCtx; 4467f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, TSMonitorVTKCtx); 4477f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(TSMonitorVTKCtx *); 4487f27e910SStefano Zampini PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKCtxCreate(const char *, TSMonitorVTKCtx *); 449fb1732b5SBarry Smith 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS); 4517cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *); 453cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec); 454e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *); 455e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *); 457d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason); 458487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *); 4595ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *); 4605ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *); 461cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *); 462cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt); 463cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *); 464cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt); 465cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool); 466dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS); 46724655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS); 4689dc50cb5SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepRollBack(TS, PetscBool *); 469d2daff3dSHong Zhang 470cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]); 471818ad0c1SBarry Smith 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal); 474e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *); 47580275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *); 47680275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal); 47780275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *); 47880275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt); 479818ad0c1SBarry Smith 48014d0ab18SJacob Faibussowitsch /*S 4818434afd1SBarry Smith TSRHSFunctionFn - A prototype of a `TS` right-hand-side evaluation function that would be passed to `TSSetRHSFunction()` 48214d0ab18SJacob Faibussowitsch 48314d0ab18SJacob Faibussowitsch Calling Sequence: 48414d0ab18SJacob Faibussowitsch + ts - timestep context 48514d0ab18SJacob Faibussowitsch . t - current time 48614d0ab18SJacob Faibussowitsch . u - input vector 48714d0ab18SJacob Faibussowitsch . F - function vector 48814d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 48914d0ab18SJacob Faibussowitsch 49014d0ab18SJacob Faibussowitsch Level: beginner 49114d0ab18SJacob Faibussowitsch 492d1c5d1fcSBarry Smith Note: 4938434afd1SBarry Smith The deprecated `TSRHSFunction` still works as a replacement for `TSRHSFunctionFn` *. 494d1c5d1fcSBarry Smith 4958434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`, 4968434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn` 49714d0ab18SJacob Faibussowitsch S*/ 4988434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSFunctionFn)(TS ts, PetscReal t, Vec u, Vec F, void *ctx); 499d1c5d1fcSBarry Smith 5008434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSFunctionFn *TSRHSFunction; 50114d0ab18SJacob Faibussowitsch 50214d0ab18SJacob Faibussowitsch /*S 5038434afd1SBarry Smith TSRHSJacobianFn - A prototype of a `TS` right-hand-side Jacobian evaluation function that would be passed to `TSSetRHSJacobian()` 50414d0ab18SJacob Faibussowitsch 50514d0ab18SJacob Faibussowitsch Calling Sequence: 50614d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 50714d0ab18SJacob Faibussowitsch . t - current time 50814d0ab18SJacob Faibussowitsch . u - input vector 50914d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix 510af27ebaaSBarry Smith . Pmat - matrix from which preconditioner is to be constructed (usually the same as `Amat`) 51114d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for matrix evaluation routine 51214d0ab18SJacob Faibussowitsch 51314d0ab18SJacob Faibussowitsch Level: beginner 51414d0ab18SJacob Faibussowitsch 515d1c5d1fcSBarry Smith Note: 5168434afd1SBarry Smith The deprecated `TSRHSJacobian` still works as a replacement for `TSRHSJacobianFn` *. 517d1c5d1fcSBarry Smith 5188434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunctionFn`, 5198434afd1SBarry Smith `TSIFunctionFn`, `TSIJacobianFn` 52014d0ab18SJacob Faibussowitsch S*/ 5218434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianFn)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx); 522d1c5d1fcSBarry Smith 5238434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianFn *TSRHSJacobian; 52414d0ab18SJacob Faibussowitsch 52514d0ab18SJacob Faibussowitsch /*S 5268434afd1SBarry Smith TSRHSJacobianPFn - A prototype of a function that computes the Jacobian of G w.r.t. the parameters P where 527af27ebaaSBarry Smith U_t = G(U,P,t), as well as the location to store the matrix that would be passed to `TSSetRHSJacobianP()` 52814d0ab18SJacob Faibussowitsch 52914d0ab18SJacob Faibussowitsch Calling Sequence: 53014d0ab18SJacob Faibussowitsch + ts - the `TS` context 53114d0ab18SJacob Faibussowitsch . t - current timestep 53214d0ab18SJacob Faibussowitsch . U - input vector (current ODE solution) 53314d0ab18SJacob Faibussowitsch . A - output matrix 53414d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 53514d0ab18SJacob Faibussowitsch 53614d0ab18SJacob Faibussowitsch Level: beginner 53714d0ab18SJacob Faibussowitsch 538d1c5d1fcSBarry Smith Note: 5398434afd1SBarry Smith The deprecated `TSRHSJacobianP` still works as a replacement for `TSRHSJacobianPFn` *. 540d1c5d1fcSBarry Smith 54114d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()` 54214d0ab18SJacob Faibussowitsch S*/ 5438434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSRHSJacobianPFn)(TS ts, PetscReal t, Vec U, Mat A, void *ctx); 54414d0ab18SJacob Faibussowitsch 5458434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSRHSJacobianPFn *TSRHSJacobianP; 546d1c5d1fcSBarry Smith 5478434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunctionFn *, void *); 5488434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunctionFn **, void **); 5498434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobianFn *, void *); 5508434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobianFn **, void **); 551e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool); 552818ad0c1SBarry Smith 55314d0ab18SJacob Faibussowitsch /*S 5548434afd1SBarry Smith TSSolutionFn - A prototype of a `TS` solution evaluation function that would be passed to `TSSetSolutionFunction()` 55514d0ab18SJacob Faibussowitsch 55614d0ab18SJacob Faibussowitsch Calling Sequence: 55714d0ab18SJacob Faibussowitsch + ts - timestep context 55814d0ab18SJacob Faibussowitsch . t - current time 55914d0ab18SJacob Faibussowitsch . u - output vector 56014d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 56114d0ab18SJacob Faibussowitsch 56214d0ab18SJacob Faibussowitsch Level: advanced 56314d0ab18SJacob Faibussowitsch 564d1c5d1fcSBarry Smith Note: 5658434afd1SBarry Smith The deprecated `TSSolutionFunction` still works as a replacement for `TSSolutionFn` *. 566d1c5d1fcSBarry Smith 56714d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()` 56814d0ab18SJacob Faibussowitsch S*/ 5698434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSSolutionFn)(TS ts, PetscReal t, Vec u, void *ctx); 57014d0ab18SJacob Faibussowitsch 5718434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSSolutionFn *TSSolutionFunction; 572d1c5d1fcSBarry Smith 5738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFn *, void *); 57414d0ab18SJacob Faibussowitsch 57514d0ab18SJacob Faibussowitsch /*S 5768434afd1SBarry Smith TSForcingFn - A prototype of a `TS` forcing function evaluation function that would be passed to `TSSetForcingFunction()` 57714d0ab18SJacob Faibussowitsch 57814d0ab18SJacob Faibussowitsch Calling Sequence: 57914d0ab18SJacob Faibussowitsch + ts - timestep context 58014d0ab18SJacob Faibussowitsch . t - current time 58114d0ab18SJacob Faibussowitsch . f - output vector 58214d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 58314d0ab18SJacob Faibussowitsch 58414d0ab18SJacob Faibussowitsch Level: advanced 58514d0ab18SJacob Faibussowitsch 586d1c5d1fcSBarry Smith Note: 5878434afd1SBarry Smith The deprecated `TSForcingFunction` still works as a replacement for `TSForcingFn` *. 588d1c5d1fcSBarry Smith 58914d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()` 59014d0ab18SJacob Faibussowitsch S*/ 5918434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSForcingFn)(TS ts, PetscReal t, Vec f, void *ctx); 59214d0ab18SJacob Faibussowitsch 5938434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSForcingFn *TSForcingFunction; 594d1c5d1fcSBarry Smith 5958434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFn *, void *); 596ef20d060SBarry Smith 59714d0ab18SJacob Faibussowitsch /*S 5988434afd1SBarry Smith TSIFunctionFn - A prototype of a `TS` implicit function evaluation function that would be passed to `TSSetIFunction() 59914d0ab18SJacob Faibussowitsch 60014d0ab18SJacob Faibussowitsch Calling Sequence: 60114d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 60214d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 60314d0ab18SJacob Faibussowitsch . U - state vector 60414d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 60514d0ab18SJacob Faibussowitsch . F - function vector 60614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function 60714d0ab18SJacob Faibussowitsch 60814d0ab18SJacob Faibussowitsch Level: beginner 60914d0ab18SJacob Faibussowitsch 610d1c5d1fcSBarry Smith Note: 6118434afd1SBarry Smith The deprecated `TSIFunction` still works as a replacement for `TSIFunctionFn` *. 612d1c5d1fcSBarry Smith 6138434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn` 61414d0ab18SJacob Faibussowitsch S*/ 6158434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIFunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec F, void *ctx); 616d1c5d1fcSBarry Smith 6178434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIFunctionFn *TSIFunction; 61814d0ab18SJacob Faibussowitsch 61914d0ab18SJacob Faibussowitsch /*S 6208434afd1SBarry Smith TSIJacobianFn - A prototype of a `TS` Jacobian evaluation function that would be passed to `TSSetIJacobian()` 62114d0ab18SJacob Faibussowitsch 62214d0ab18SJacob Faibussowitsch Calling Sequence: 62314d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 62414d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 62514d0ab18SJacob Faibussowitsch . U - state vector 62614d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 62714d0ab18SJacob Faibussowitsch . a - shift 62814d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 62914d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat` 63014d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for Jacobian evaluation routine 63114d0ab18SJacob Faibussowitsch 63214d0ab18SJacob Faibussowitsch Level: beginner 63314d0ab18SJacob Faibussowitsch 634d1c5d1fcSBarry Smith Note: 6358434afd1SBarry Smith The deprecated `TSIJacobian` still works as a replacement for `TSIJacobianFn` *. 63614d0ab18SJacob Faibussowitsch 6378434afd1SBarry Smith .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunctionFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn` 638d1c5d1fcSBarry Smith S*/ 6398434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSIJacobianFn)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx); 640d1c5d1fcSBarry Smith 6418434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSIJacobianFn *TSIJacobian; 642d1c5d1fcSBarry Smith 6438434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunctionFn *, void *); 6448434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunctionFn **, void **); 6458434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobianFn *, void *); 6468434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobianFn **, void **); 647316643e7SJed Brown 64814d0ab18SJacob Faibussowitsch /*S 6498434afd1SBarry Smith TSI2FunctionFn - A prototype of a `TS` implicit function evaluation function for 2nd order systems that would be passed to `TSSetI2Function()` 65014d0ab18SJacob Faibussowitsch 65114d0ab18SJacob Faibussowitsch Calling Sequence: 65214d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 65314d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 65414d0ab18SJacob Faibussowitsch . U - state vector 65514d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 65614d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector 65714d0ab18SJacob Faibussowitsch . F - function vector 65814d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for matrix evaluation routine (may be `NULL`) 65914d0ab18SJacob Faibussowitsch 66014d0ab18SJacob Faibussowitsch Level: advanced 66114d0ab18SJacob Faibussowitsch 662d1c5d1fcSBarry Smith Note: 6638434afd1SBarry Smith The deprecated `TSI2Function` still works as a replacement for `TSI2FunctionFn` *. 664d1c5d1fcSBarry Smith 6658434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()`, `TSIFunctionFn` 66614d0ab18SJacob Faibussowitsch S*/ 6678434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSI2FunctionFn)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx); 668d1c5d1fcSBarry Smith 6698434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2FunctionFn *TSI2Function; 67014d0ab18SJacob Faibussowitsch 67114d0ab18SJacob Faibussowitsch /*S 6728434afd1SBarry Smith TSI2JacobianFn - A prototype of a `TS` implicit Jacobian evaluation function for 2nd order systems that would be passed to `TSSetI2Jacobian()` 67314d0ab18SJacob Faibussowitsch 67414d0ab18SJacob Faibussowitsch Calling Sequence: 67514d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 67614d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 67714d0ab18SJacob Faibussowitsch . U - state vector 67814d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 67914d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector 68014d0ab18SJacob Faibussowitsch . v - shift for U_t 68114d0ab18SJacob Faibussowitsch . a - shift for U_tt 68214d0ab18SJacob Faibussowitsch . 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 68314d0ab18SJacob Faibussowitsch . jac - matrix from which to construct the preconditioner, may be same as `J` 68414d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for matrix evaluation routine 68514d0ab18SJacob Faibussowitsch 68614d0ab18SJacob Faibussowitsch Level: advanced 68714d0ab18SJacob Faibussowitsch 688d1c5d1fcSBarry Smith Note: 6898434afd1SBarry Smith The deprecated `TSI2Jacobian` still works as a replacement for `TSI2JacobianFn` *. 69014d0ab18SJacob Faibussowitsch 6918434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunctionFn`, `TSIJacobianFn`, `TSRHSFunctionFn`, `TSRHSJacobianFn` 692d1c5d1fcSBarry Smith S*/ 6938434afd1SBarry Smith 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); 694d1c5d1fcSBarry Smith 6958434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSI2JacobianFn *TSI2Jacobian; 696d1c5d1fcSBarry Smith 6978434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2FunctionFn *, void *); 6988434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2FunctionFn **, void **); 6998434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2JacobianFn *, void *); 7008434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2JacobianFn **, void **); 701efe9872eSLisandro Dalcin 702545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS); 703545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *); 7048434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunctionFn *, void *); 7053a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIFunction(TS, const char[], Vec, TSIFunctionFn *, void *); 7063a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIJacobian(TS, const char[], Mat, Mat, TSIJacobianFn *, void *); 7071d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *); 7081d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]); 7090fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool); 7100fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *); 7114bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSNES(TS, SNES *); 7124bd3aaa3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetSNES(TS, SNES); 713d9194312SHong Zhang 7148434afd1SBarry Smith PETSC_EXTERN TSRHSFunctionFn TSComputeRHSFunctionLinear; 7158434afd1SBarry Smith PETSC_EXTERN TSRHSJacobianFn TSComputeRHSJacobianConstant; 716014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *); 717d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 7181c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec); 7199b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec); 720847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 721d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat); 722e34be4c2SBarry Smith 723014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 724b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal)); 7259be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *)); 726c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS)); 727014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 728ecc87898SStefano Zampini PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscBool, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *); 729014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS); 730b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal); 7319be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *); 732c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS); 733014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS); 7346bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS); 7356bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *); 7366bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec); 737c6bf8827SStefano Zampini PETSC_EXTERN PetscErrorCode TSGetStepResize(TS, PetscBool *); 7386bd3a4fdSStefano Zampini 739014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec); 740014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec); 741c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *); 7427453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *); 7438a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *); 744014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal); 745014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *); 746d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *)); 747d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *); 748000e7ae3SMatthew Knepley 749014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *); 7508d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *); 751014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *); 752014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal); 753014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *); 7548d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *); 755014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *); 756014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal); 757014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 75821c89e3eSBarry Smith 759014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]); 760ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]); 7611d6018f0SLisandro Dalcin 762014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec); 763d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat); 764014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool); 765d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool); 766efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec); 767efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat); 768f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *); 769818ad0c1SBarry Smith 770014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec); 771d6ebe24aSShri Abhyankar 772967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 7738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunctionFn *, void *); 7748434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunctionFn **, void **); 77549abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscCtxDestroyFn *); 7768434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobianFn *, void *); 7778434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobianFn **, void **); 77849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscCtxDestroyFn *); 7798434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunctionFn *, void *); 7808434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunctionFn **, void **); 78149abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscCtxDestroyFn *); 7828434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobianFn *, void *); 7838434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobianFn **, void **); 78449abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscCtxDestroyFn *); 7858434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2FunctionFn *, void *); 7868434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2FunctionFn **, void **); 78749abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscCtxDestroyFn *); 7888434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2JacobianFn *, void *); 7898434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2JacobianFn **, void **); 79049abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscCtxDestroyFn *); 791efe9872eSLisandro Dalcin 79214d0ab18SJacob Faibussowitsch /*S 7938434afd1SBarry Smith TSTransientVariableFn - A prototype of a function to transform from state to transient variables that would be passed to `TSSetTransientVariable()` 79414d0ab18SJacob Faibussowitsch 79514d0ab18SJacob Faibussowitsch Calling Sequence: 79614d0ab18SJacob Faibussowitsch + ts - timestep context 79714d0ab18SJacob Faibussowitsch . p - input vector (primitive form) 79814d0ab18SJacob Faibussowitsch . c - output vector, transient variables (conservative form) 79914d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 80014d0ab18SJacob Faibussowitsch 80114d0ab18SJacob Faibussowitsch Level: advanced 80214d0ab18SJacob Faibussowitsch 803d1c5d1fcSBarry Smith Note: 8048434afd1SBarry Smith The deprecated `TSTransientVariable` still works as a replacement for `TSTransientVariableFn` *. 805d1c5d1fcSBarry Smith 80614d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()` 80714d0ab18SJacob Faibussowitsch S*/ 8088434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSTransientVariableFn)(TS ts, Vec p, Vec c, void *ctx); 80914d0ab18SJacob Faibussowitsch 8108434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSTransientVariableFn *TSTransientVariable; 811d1c5d1fcSBarry Smith 8128434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariableFn *, void *); 8138434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariableFn *, void *); 8148434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariableFn **, void *); 815e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec); 816e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *); 817e3c11fc1SJed Brown 8188434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFn *, void *); 8198434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFn **, void **); 8208434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFn *, void *); 8218434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFn **, void **); 822e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *); 823e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *); 8247f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec); 8259b7cd975SBarry Smith 8260c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **); 827d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *); 8280c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **); 829d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *); 8300c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **); 831d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 832b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM); 833b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM); 834b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM); 8356c6b9e74SPeter Brune 8366c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer)); 8376c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer)); 8386c6b9e74SPeter Brune 83914d0ab18SJacob Faibussowitsch /*S 840dd8e379bSPierre Jolivet DMDATSRHSFunctionLocalFn - A prototype of a local `TS` right-hand side residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSFunctionLocal()` 8416c6b9e74SPeter Brune 84214d0ab18SJacob Faibussowitsch Calling Sequence: 84314d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 84414d0ab18SJacob Faibussowitsch . t - time at which to evaluate residual 84514d0ab18SJacob Faibussowitsch . x - array of local state information 84614d0ab18SJacob Faibussowitsch . f - output array of local residual information 84714d0ab18SJacob Faibussowitsch - ctx - optional user context 84814d0ab18SJacob Faibussowitsch 84914d0ab18SJacob Faibussowitsch Level: beginner 85014d0ab18SJacob Faibussowitsch 851d1c5d1fcSBarry Smith Note: 8528434afd1SBarry Smith The deprecated `DMDATSRHSFunctionLocal` still works as a replacement for `DMDATSRHSFunctionLocalFn` *. 853d1c5d1fcSBarry Smith 8548434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunctionFn`, `DMDATSRHSJacobianLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn` 85514d0ab18SJacob Faibussowitsch S*/ 8568434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx); 857d1c5d1fcSBarry Smith 8588434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSFunctionLocalFn *DMDATSRHSFunctionLocal; 85914d0ab18SJacob Faibussowitsch 86014d0ab18SJacob Faibussowitsch /*S 8618434afd1SBarry Smith DMDATSRHSJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetRHSJacobianLocal()` 86214d0ab18SJacob Faibussowitsch 86314d0ab18SJacob Faibussowitsch Calling Sequence: 86414d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 86514d0ab18SJacob Faibussowitsch . t - time at which to evaluate residual 86614d0ab18SJacob Faibussowitsch . x - array of local state information 86714d0ab18SJacob Faibussowitsch . J - Jacobian matrix 86814d0ab18SJacob Faibussowitsch . B - matrix from which to construct the preconditioner; often same as `J` 86914d0ab18SJacob Faibussowitsch - ctx - optional context 87014d0ab18SJacob Faibussowitsch 87114d0ab18SJacob Faibussowitsch Level: beginner 87214d0ab18SJacob Faibussowitsch 873d1c5d1fcSBarry Smith Note: 8748434afd1SBarry Smith The deprecated `DMDATSRHSJacobianLocal` still works as a replacement for `DMDATSRHSJacobianLocalFn` *. 875d1c5d1fcSBarry Smith 8768434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobianFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSIJacobianLocalFn`, `DMDATSIFunctionLocalFn` 87714d0ab18SJacob Faibussowitsch S*/ 8788434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSRHSJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx); 879d1c5d1fcSBarry Smith 8808434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSRHSJacobianLocalFn *DMDATSRHSJacobianLocal; 88114d0ab18SJacob Faibussowitsch 88214d0ab18SJacob Faibussowitsch /*S 8838434afd1SBarry Smith DMDATSIFunctionLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIFunctionLocal()` 88414d0ab18SJacob Faibussowitsch 88514d0ab18SJacob Faibussowitsch Calling Sequence: 88614d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 88714d0ab18SJacob Faibussowitsch . t - time at which to evaluate residual 88814d0ab18SJacob Faibussowitsch . x - array of local state information 88914d0ab18SJacob Faibussowitsch . xdot - array of local time derivative information 89014d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information 89114d0ab18SJacob Faibussowitsch - ctx - optional context 89214d0ab18SJacob Faibussowitsch 89314d0ab18SJacob Faibussowitsch Level: beginner 89414d0ab18SJacob Faibussowitsch 895d1c5d1fcSBarry Smith Note: 8968434afd1SBarry Smith The deprecated `DMDATSIFunctionLocal` still works as a replacement for `DMDATSIFunctionLocalFn` *. 897d1c5d1fcSBarry Smith 8988434afd1SBarry Smith .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocalFn`, `TSIFunctionFn` 89914d0ab18SJacob Faibussowitsch S*/ 9008434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIFunctionLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx); 901d1c5d1fcSBarry Smith 9028434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIFunctionLocalFn *DMDATSIFunctionLocal; 90314d0ab18SJacob Faibussowitsch 90414d0ab18SJacob Faibussowitsch /*S 9058434afd1SBarry Smith DMDATSIJacobianLocalFn - A prototype of a local residual evaluation function for use with `DMDA` that would be passed to `DMDATSSetIJacobianLocal()` 90614d0ab18SJacob Faibussowitsch 90714d0ab18SJacob Faibussowitsch Calling Sequence: 90814d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 90914d0ab18SJacob Faibussowitsch . t - time at which to evaluate the jacobian 91014d0ab18SJacob Faibussowitsch . x - array of local state information 91114d0ab18SJacob Faibussowitsch . xdot - time derivative at this state 91214d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter 91314d0ab18SJacob Faibussowitsch . J - Jacobian matrix 91414d0ab18SJacob Faibussowitsch . B - matrix from which to construct the preconditioner; often same as `J` 91514d0ab18SJacob Faibussowitsch - ctx - optional context 91614d0ab18SJacob Faibussowitsch 91714d0ab18SJacob Faibussowitsch Level: beginner 91814d0ab18SJacob Faibussowitsch 919d1c5d1fcSBarry Smith Note: 9208434afd1SBarry Smith The deprecated `DMDATSIJacobianLocal` still works as a replacement for `DMDATSIJacobianLocalFn` *. 92114d0ab18SJacob Faibussowitsch 9228434afd1SBarry Smith .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobianFn`, `DMDATSIFunctionLocalFn`, `DMDATSRHSFunctionLocalFn`, `DMDATSRHSJacobianlocal()` 923d1c5d1fcSBarry Smith S*/ 9248434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(DMDATSIJacobianLocalFn)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx); 925d1c5d1fcSBarry Smith 9268434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef DMDATSIJacobianLocalFn *DMDATSIJacobianLocal; 927d1c5d1fcSBarry Smith 9288434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocalFn *, void *); 9298434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocalFn *, void *); 9308434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocalFn *, void *); 9318434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocalFn *, void *); 9326c6b9e74SPeter Brune 933be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx; 934d1212d36SBarry Smith typedef struct { 935d1212d36SBarry Smith Vec ray; 936d1212d36SBarry Smith VecScatter scatter; 937d1212d36SBarry Smith PetscViewer viewer; 93851b4a12fSMatthew G. Knepley TSMonitorLGCtx lgctx; 939d1212d36SBarry Smith } TSMonitorDMDARayCtx; 940d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **); 941d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *); 94251b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *); 943d1212d36SBarry Smith 944bdad233fSMatthew Knepley /* Dynamic creation and loading functions */ 945140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList; 94619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSGetType(TS, TSType *); 94719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSSetType(TS, TSType); 948bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 94930de9b25SBarry Smith 950014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *); 951deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES); 952014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *); 953818ad0c1SBarry Smith 954014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer); 95555849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer); 956fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]); 957fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]); 95855849f57SBarry Smith 95955849f57SBarry Smith #define TS_FILE_CLASSID 1211225 96021c89e3eSBarry Smith 961014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *); 962014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *); 96321c89e3eSBarry Smith 964a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *); 965a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *); 9664f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *); 9674f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *); 96831152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *); 969b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **); 970e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *); 971a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *); 972a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *); 97349abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *); 97449abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscCtxDestroyFn *, void *); 9754f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *); 976201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *); 977201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *); 978edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 979d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 980ef20d060SBarry Smith 981e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork { 982e669de00SBarry Smith PetscInt nlg; 983e669de00SBarry Smith PetscDrawLG *lg; 984e669de00SBarry Smith PetscBool semilogy; 985e669de00SBarry Smith PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */ 986e669de00SBarry Smith }; 987e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork; 988e669de00SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *); 989e669de00SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *); 990e669de00SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *); 991e669de00SBarry Smith 992b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx; 993b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *); 994b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *); 995b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *); 996b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *); 997b3d3934dSBarry Smith 9988189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx; 9998189c53fSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *); 10008189c53fSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *); 10018189c53fSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *); 10023914022bSBarry Smith 10031b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx; 100460e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *); 10051b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *); 10060ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 10071b575b74SJoseph Pusztay 100860e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx; 100960e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *); 101060e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 101160e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *); 101260e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 101360e16b1bSMatthew G. Knepley 1014ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *); 1015ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal); 1016fe4ad979SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal); 1017fe4ad979SIlya Fursov PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt) 1018fe4ad979SIlya Fursov { 1019fe4ad979SIlya Fursov return TSSetPostEventSecondStep(ts, dt); 1020fe4ad979SIlya Fursov } 10216427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]); 10221ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *); 10231ea83e56SMiguel 102476bdecfbSBarry Smith /*J 1025195e9b02SBarry Smith TSSSPType - string with the name of a `TSSSP` scheme. 1026815f1ad5SJed Brown 1027815f1ad5SJed Brown Level: beginner 1028815f1ad5SJed Brown 10291cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP` 103076bdecfbSBarry Smith J*/ 103119fd82e9SBarry Smith typedef const char *TSSSPType; 1032815f1ad5SJed Brown #define TSSSPRKS2 "rks2" 1033815f1ad5SJed Brown #define TSSSPRKS3 "rks3" 1034815f1ad5SJed Brown #define TSSSPRK104 "rk104" 1035815f1ad5SJed Brown 103619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSSSPSetType(TS, TSSSPType); 103719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSSSPGetType(TS, TSSSPType *); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS, PetscInt); 1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS, PetscInt *); 1040787849ffSJed Brown PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 1041560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 1042013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList; 1043815f1ad5SJed Brown 1044ed657a08SJed Brown /*S 104584df9cb4SJed Brown TSAdapt - Abstract object that manages time-step adaptivity 104684df9cb4SJed Brown 104784df9cb4SJed Brown Level: beginner 104884df9cb4SJed Brown 1049af27ebaaSBarry Smith .seealso: [](ch_ts), `TS`, `TSGetAdapt()`, `TSAdaptCreate()`, `TSAdaptType` 105084df9cb4SJed Brown S*/ 105184df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt; 105284df9cb4SJed Brown 1053f0fc11ceSJed Brown /*J 105487497f52SBarry Smith TSAdaptType - String with the name of `TSAdapt` scheme. 105584df9cb4SJed Brown 105684df9cb4SJed Brown Level: beginner 105784df9cb4SJed Brown 1058af27ebaaSBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptSetType()`, `TS`, `TSAdapt` 1059f0fc11ceSJed Brown J*/ 1060f8963224SJed Brown typedef const char *TSAdaptType; 1061b0f836d7SJed Brown #define TSADAPTNONE "none" 10621917a363SLisandro Dalcin #define TSADAPTBASIC "basic" 1063cb7ab186SLisandro Dalcin #define TSADAPTDSP "dsp" 10648d59e960SJed Brown #define TSADAPTCFL "cfl" 10651917a363SLisandro Dalcin #define TSADAPTGLEE "glee" 1066d949e4a4SStefano Zampini #define TSADAPTHISTORY "history" 106784df9cb4SJed Brown 1068552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *); 1069bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt)); 1070607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *); 107319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType); 1074d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *); 1075014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]); 1076014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 1077014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool); 1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **); 1079014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *); 1080b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *); 1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer); 1082ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer); 1083dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *); 1084099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt); 1085014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *); 1086014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool); 1087bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool); 10881917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal); 10891917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *); 109076cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal); 109176cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *); 10921917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal); 10931917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *); 109462c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal); 109562c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *); 1096014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal); 10971917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *); 1098b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *)); 1099d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool); 110075017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool); 110175017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *); 1102de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt); 1103cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *); 1104cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal); 110584df9cb4SJed Brown 110684df9cb4SJed Brown /*S 110787497f52SBarry Smith TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE` 1108ed657a08SJed Brown 1109ed657a08SJed Brown Level: beginner 1110ed657a08SJed Brown 1111195e9b02SBarry Smith Developer Note: 111287497f52SBarry Smith This functionality should be replaced by the `TSAdapt`. 111384df9cb4SJed Brown 11141cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType` 1115ed657a08SJed Brown S*/ 111626d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt; 1117ed657a08SJed Brown 111876bdecfbSBarry Smith /*J 111987497f52SBarry Smith TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme 1120ed657a08SJed Brown 1121ed657a08SJed Brown Level: beginner 1122ed657a08SJed Brown 1123195e9b02SBarry Smith Developer Note: 1124195e9b02SBarry Smith This functionality should be replaced by the `TSAdaptType`. 1125195e9b02SBarry Smith 11261cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS` 112776bdecfbSBarry Smith J*/ 112826d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType; 112926d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none" 113026d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size" 113126d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both" 1132ed657a08SJed Brown 113326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt)); 113426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void); 113526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void); 113626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *); 113726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType); 113826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]); 113926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *); 114026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer); 1141dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *); 114226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *); 1143ed657a08SJed Brown 114476bdecfbSBarry Smith /*J 114587497f52SBarry Smith TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme 1146ed657a08SJed Brown 1147ed657a08SJed Brown Level: beginner 1148ed657a08SJed Brown 11491cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept` 115076bdecfbSBarry Smith J*/ 115126d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType; 115226d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always" 1153ed657a08SJed Brown 1154d1c5d1fcSBarry Smith /*S 11558434afd1SBarry Smith TSGLLEAcceptFn - A prototype of a `TS` accept function that would be passed to `TSGLLEAcceptRegister()` 1156d1c5d1fcSBarry Smith 1157d1c5d1fcSBarry Smith Calling Sequence: 1158d1c5d1fcSBarry Smith + ts - timestep context 1159d1c5d1fcSBarry Smith . nt - time to end of solution time 1160d1c5d1fcSBarry Smith . h - the proposed step-size 1161d1c5d1fcSBarry Smith . enorm - unknown 1162d1c5d1fcSBarry Smith - accept - output, if the proposal is accepted 1163d1c5d1fcSBarry Smith 1164d1c5d1fcSBarry Smith Level: beginner 1165d1c5d1fcSBarry Smith 1166d1c5d1fcSBarry Smith Note: 11678434afd1SBarry Smith The deprecated `TSGLLEAcceptFunction` still works as a replacement for `TSGLLEAcceptFn` * 1168d1c5d1fcSBarry Smith 11698434afd1SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunctionFn`, 11708434afd1SBarry Smith `TSIJacobianFn`, `TSRHSJacobianFn`, `TSGLLEAcceptRegister()` 1171d1c5d1fcSBarry Smith S*/ 11728434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSGLLEAcceptFn)(TS ts, PetscReal nt, PetscReal h, const PetscReal enorm[], PetscBool *accept); 1173d1c5d1fcSBarry Smith 11748434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSGLLEAcceptFn *TSGLLEAcceptFunction; 1175d1c5d1fcSBarry Smith 11768434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFn *); 1177ed657a08SJed Brown 117876bdecfbSBarry Smith /*J 1179195e9b02SBarry Smith TSGLLEType - string with the name of a General Linear `TSGLLE` type 118018b56cb9SJed Brown 118118b56cb9SJed Brown Level: beginner 118218b56cb9SJed Brown 11831cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept` 118476bdecfbSBarry Smith J*/ 118526d28e4eSEmil Constantinescu typedef const char *TSGLLEType; 118626d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks" 118718b56cb9SJed Brown 118826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS)); 118926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void); 119026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void); 119126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType); 119226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *); 119326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType); 119418b56cb9SJed Brown 1195020d8f30SJed Brown /*J 1196195e9b02SBarry Smith TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type 11977d5697caSHong Zhang 11987d5697caSHong Zhang Level: beginner 11997d5697caSHong Zhang 12001cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()` 12017d5697caSHong Zhang J*/ 12027d5697caSHong Zhang #define TSEIMEXType char * 12037d5697caSHong Zhang 12047d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt); 12057d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt); 12067d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool); 12077d5697caSHong Zhang 12087d5697caSHong Zhang /*J 1209195e9b02SBarry Smith TSRKType - String with the name of a Runge-Kutta `TSRK` type 1210f68a32c8SEmil Constantinescu 1211f68a32c8SEmil Constantinescu Level: beginner 1212f68a32c8SEmil Constantinescu 12131cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()` 1214f68a32c8SEmil Constantinescu J*/ 1215f68a32c8SEmil Constantinescu typedef const char *TSRKType; 1216f68a32c8SEmil Constantinescu #define TSRK1FE "1fe" 1217fdefd5e5SDebojyoti Ghosh #define TSRK2A "2a" 1218e7685601SHong Zhang #define TSRK2B "2b" 1219f68a32c8SEmil Constantinescu #define TSRK3 "3" 1220fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs" 1221f68a32c8SEmil Constantinescu #define TSRK4 "4" 1222fdefd5e5SDebojyoti Ghosh #define TSRK5F "5f" 1223fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp" 122405e23783SLisandro Dalcin #define TSRK5BS "5bs" 122537acaa02SHendrik Ranocha #define TSRK6VR "6vr" 122637acaa02SHendrik Ranocha #define TSRK7VR "7vr" 122737acaa02SHendrik Ranocha #define TSRK8VR "8vr" 122805e23783SLisandro Dalcin 12292ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *); 12300f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *); 12310f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType); 12320f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *); 12330fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool); 12340fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *); 1235f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1236f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 1237560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 1238f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 1239f68a32c8SEmil Constantinescu 1240f68a32c8SEmil Constantinescu /*J 1241af27ebaaSBarry Smith TSMPRKType - String with the name of a partitioned Runge-Kutta `TSMPRK` type 12429f537349Sluzhanghpp 12439f537349Sluzhanghpp Level: beginner 12449f537349Sluzhanghpp 12451cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()` 12469f537349Sluzhanghpp J*/ 12474b84eec9SHong Zhang typedef const char *TSMPRKType; 124819c2959aSHong Zhang #define TSMPRK2A22 "2a22" 124919c2959aSHong Zhang #define TSMPRK2A23 "2a23" 125019c2959aSHong Zhang #define TSMPRK2A32 "2a32" 125119c2959aSHong Zhang #define TSMPRK2A33 "2a33" 12524b84eec9SHong Zhang #define TSMPRKP2 "p2" 12534b84eec9SHong Zhang #define TSMPRKP3 "p3" 12549f537349Sluzhanghpp 12554b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *); 12564b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType); 125779bc8a77SHong Zhang 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[]); 12584b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void); 12594b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void); 12604b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void); 12619f537349Sluzhanghpp 12629f537349Sluzhanghpp /*J 1263195e9b02SBarry Smith TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type 1264d2567f34SHong Zhang 1265d2567f34SHong Zhang Level: beginner 1266d2567f34SHong Zhang 12671cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()` 1268d2567f34SHong Zhang J*/ 1269d2567f34SHong Zhang typedef const char *TSIRKType; 1270d2567f34SHong Zhang #define TSIRKGAUSS "gauss" 1271d2567f34SHong Zhang 1272d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *); 1273d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType); 1274d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *); 1275d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt); 1276d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS)); 1277d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *); 1278d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void); 1279d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void); 1280d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void); 1281d2567f34SHong Zhang 1282d2567f34SHong Zhang /*J 1283195e9b02SBarry Smith TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type 1284b6a60446SDebojyoti Ghosh 1285b6a60446SDebojyoti Ghosh Level: beginner 1286b6a60446SDebojyoti Ghosh 12871cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()` 1288b6a60446SDebojyoti Ghosh J*/ 1289b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType; 1290ab8c5912SEmil Constantinescu #define TSGLEEi1 "BE1" 1291b6a60446SDebojyoti Ghosh #define TSGLEE23 "23" 1292b6a60446SDebojyoti Ghosh #define TSGLEE24 "24" 1293b6a60446SDebojyoti Ghosh #define TSGLEE25I "25i" 1294b6a60446SDebojyoti Ghosh #define TSGLEE35 "35" 1295b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A "exrk2a" 1296b6a60446SDebojyoti Ghosh #define TSGLEERK32G1 "rk32g1" 1297b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex" 129816a05f60SBarry Smith 1299b6a60446SDebojyoti Ghosh /*J 1300195e9b02SBarry Smith TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type 1301b6a60446SDebojyoti Ghosh 1302b6a60446SDebojyoti Ghosh Level: beginner 1303b6a60446SDebojyoti Ghosh 13041cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()` 1305b6a60446SDebojyoti Ghosh J*/ 1306b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *); 1307b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType); 130857df6a1bSDebojyoti Ghosh 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[]); 13094bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void); 1310b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void); 1311b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void); 1312b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void); 1313b6a60446SDebojyoti Ghosh 1314b6a60446SDebojyoti Ghosh /*J 1315195e9b02SBarry Smith TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type 1316020d8f30SJed Brown 1317020d8f30SJed Brown Level: beginner 1318020d8f30SJed Brown 13191cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()` 1320020d8f30SJed Brown J*/ 132119fd82e9SBarry Smith typedef const char *TSARKIMEXType; 1322e817cc15SEmil Constantinescu #define TSARKIMEX1BEE "1bee" 13231f80e275SEmil Constantinescu #define TSARKIMEXA2 "a2" 13241f80e275SEmil Constantinescu #define TSARKIMEXL2 "l2" 13251f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122" 13261f80e275SEmil Constantinescu #define TSARKIMEX2C "2c" 13278a381b04SJed Brown #define TSARKIMEX2D "2d" 1328a3a57f36SJed Brown #define TSARKIMEX2E "2e" 13296cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2" 13308a381b04SJed Brown #define TSARKIMEX3 "3" 13316cf0794eSJed Brown #define TSARKIMEXBPR3 "bpr3" 13326cf0794eSJed Brown #define TSARKIMEXARS443 "ars443" 13338a381b04SJed Brown #define TSARKIMEX4 "4" 13348a381b04SJed Brown #define TSARKIMEX5 "5" 133519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *); 133619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType); 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool); 13383a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *); 13393a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXSetFastSlowSplit(TS, PetscBool); 13403a2a065fSHong Zhang PETSC_EXTERN PetscErrorCode TSARKIMEXGetFastSlowSplit(TS, PetscBool *); 134119fd82e9SBarry Smith 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[]); 1342607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 1343560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 1344014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 13458a381b04SJed Brown 1346020d8f30SJed Brown /*J 1347d27334e2SStefano Zampini TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type 1348d27334e2SStefano Zampini 1349d27334e2SStefano Zampini Level: beginner 1350d27334e2SStefano Zampini 1351d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()` 1352d27334e2SStefano Zampini J*/ 1353d27334e2SStefano Zampini typedef const char *TSDIRKType; 1354d27334e2SStefano Zampini #define TSDIRKS212 "s212" 13553405e92cSStefano Zampini #define TSDIRKES122SAL "es122sal" 1356d27334e2SStefano Zampini #define TSDIRKES213SAL "es213sal" 1357d27334e2SStefano Zampini #define TSDIRKES324SAL "es324sal" 1358d27334e2SStefano Zampini #define TSDIRKES325SAL "es325sal" 1359d27334e2SStefano Zampini #define TSDIRK657A "657a" 1360d27334e2SStefano Zampini #define TSDIRKES648SA "es648sa" 1361d27334e2SStefano Zampini #define TSDIRK658A "658a" 1362d27334e2SStefano Zampini #define TSDIRKS659A "s659a" 1363d27334e2SStefano Zampini #define TSDIRK7510SAL "7510sal" 1364d27334e2SStefano Zampini #define TSDIRKES7510SA "es7510sa" 1365d27334e2SStefano Zampini #define TSDIRK759A "759a" 1366d27334e2SStefano Zampini #define TSDIRKS7511SAL "s7511sal" 1367d27334e2SStefano Zampini #define TSDIRK8614A "8614a" 1368d27334e2SStefano Zampini #define TSDIRK8616SAL "8616sal" 1369d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal" 1370d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *); 1371d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType); 1372d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1373d27334e2SStefano Zampini 1374d27334e2SStefano Zampini /*J 1375195e9b02SBarry Smith TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type 1376020d8f30SJed Brown 1377020d8f30SJed Brown Level: beginner 1378020d8f30SJed Brown 13791cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()` 1380020d8f30SJed Brown J*/ 138119fd82e9SBarry Smith typedef const char *TSRosWType; 138261692a83SJed Brown #define TSROSW2M "2m" 138361692a83SJed Brown #define TSROSW2P "2p" 1384fe7e6d57SJed Brown #define TSROSWRA3PW "ra3pw" 1385fe7e6d57SJed Brown #define TSROSWRA34PW2 "ra34pw2" 1386337ef93bSJed Brown #define TSROSWR34PRW "r34prw" 1387337ef93bSJed Brown #define TSROSWR3PRL2 "r3prl2" 1388ef3c5b88SJed Brown #define TSROSWRODAS3 "rodas3" 138948665b53SJed Brown #define TSROSWRODASPR "rodaspr" 1390337ef93bSJed Brown #define TSROSWRODASPR2 "rodaspr2" 1391ef3c5b88SJed Brown #define TSROSWSANDU3 "sandu3" 1392961f28d0SJed Brown #define TSROSWASSP3P3S1C "assp3p3s1c" 1393961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c" 139443b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 1395753f8adbSEmil Constantinescu #define TSROSWARK3 "ark3" 13963606a31eSEmil Constantinescu #define TSROSWTHETA1 "theta1" 13973606a31eSEmil Constantinescu #define TSROSWTHETA2 "theta2" 139842faf41dSJed Brown #define TSROSWGRK4T "grk4t" 139942faf41dSJed Brown #define TSROSWSHAMP4 "shamp4" 140042faf41dSJed Brown #define TSROSWVELDD4 "veldd4" 140142faf41dSJed Brown #define TSROSW4L "4l" 1402b1c69cc3SEmil Constantinescu 14036bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *); 14046bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType); 1405014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool); 140619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 140719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal); 1408607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 1409560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 1410014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 1411e27a552bSJed Brown 1412211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt); 1413211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *); 1414211a84d6SLisandro Dalcin 14156bd7aeb5SHong Zhang /*J 1416195e9b02SBarry Smith TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type 14176bd7aeb5SHong Zhang 14186bd7aeb5SHong Zhang Level: beginner 14196bd7aeb5SHong Zhang 14201cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()` 14216bd7aeb5SHong Zhang J*/ 1422e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType; 1423e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER "1" 1424e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2" 1425e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3 "3" 1426e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4 "4" 1427e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType); 1428e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *); 1429e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]); 14304bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void); 1431e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void); 1432e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void); 1433e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void); 14346bd7aeb5SHong Zhang 1435db4653c2SDaniel Finn /*J 1436195e9b02SBarry Smith TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy), 143787497f52SBarry Smith but also has the property for some systems of monotonicity in a functional. 1438db4653c2SDaniel Finn 1439db4653c2SDaniel Finn Level: beginner 1440db4653c2SDaniel Finn 14411cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()` 1442db4653c2SDaniel Finn J*/ 144340be0ff1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSDiscGradSetFormulation(TS, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*)(TS, PetscReal, Vec, Vec, void *), void *); 14444bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSDiscGradGetFormulation(TS, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (**)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (**)(TS, PetscReal, Vec, Vec, void *), void *); 1445db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *); 1446db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool); 144740be0ff1SMatthew G. Knepley 144883e2fdc7SBarry Smith /* 14490f3b3ca1SHong Zhang PETSc interface to Sundials 145083e2fdc7SBarry Smith */ 1451e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2 14529371c9d4SSatish Balay typedef enum { 14539371c9d4SSatish Balay SUNDIALS_ADAMS = 1, 14549371c9d4SSatish Balay SUNDIALS_BDF = 2 14559371c9d4SSatish Balay } TSSundialsLmmType; 14566a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 14579371c9d4SSatish Balay typedef enum { 14589371c9d4SSatish Balay SUNDIALS_MODIFIED_GS = 1, 14599371c9d4SSatish Balay SUNDIALS_CLASSICAL_GS = 2 14609371c9d4SSatish Balay } TSSundialsGramSchmidtType; 14616a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 14624bf303faSJacob Faibussowitsch 1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType); 1464014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *); 1465014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal); 1466014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal); 1467014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal); 1468014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *); 1469014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType); 1470014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt); 1471014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal); 1472014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool); 1473014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt); 1474c4e80e11SFlorian PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt); 1475918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool); 14766dc17ff5SMatthew Knepley #endif 147783e2fdc7SBarry Smith 1478014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal); 1479014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *); 1480014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *); 1481014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool); 14820de4c49aSJed Brown 1483014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal); 1484014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal); 1485014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *); 148688df8a41SLisandro Dalcin 1487818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal); 1488818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal); 1489818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 1490818efac9SLisandro Dalcin 1491220f924aSDavid Kamensky /*S 14928434afd1SBarry Smith TSAlpha2PredictorFn - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in 1493220f924aSDavid Kamensky a second-order generalized-alpha time integrator. 1494220f924aSDavid Kamensky 1495220f924aSDavid Kamensky Calling Sequence: 1496220f924aSDavid Kamensky + ts - the `TS` context obtained from `TSCreate()` 1497220f924aSDavid Kamensky . X0 - the previous time step's state vector 1498220f924aSDavid Kamensky . V0 - the previous time step's first derivative of the state vector 1499220f924aSDavid Kamensky . A0 - the previous time step's second derivative of the state vector 1500220f924aSDavid Kamensky . X1 - the vector into which the initial guess for the current time step will be written 1501220f924aSDavid Kamensky - ctx - [optional] user-defined context for the predictor evaluation routine (may be `NULL`) 1502220f924aSDavid Kamensky 1503220f924aSDavid Kamensky Level: intermediate 1504220f924aSDavid Kamensky 1505d1c5d1fcSBarry Smith Note: 15068434afd1SBarry Smith The deprecated `TSAlpha2Predictor` still works as a replacement for `TSAlpha2PredictorFn` *. 1507d1c5d1fcSBarry Smith 1508220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()` 1509220f924aSDavid Kamensky S*/ 15108434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(TSAlpha2PredictorFn)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx); 1511d1c5d1fcSBarry Smith 15128434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef TSAlpha2PredictorFn *TSAlpha2Predictor; 1513d1c5d1fcSBarry Smith 15148434afd1SBarry Smith PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2PredictorFn *, void *ctx); 1515220f924aSDavid Kamensky 1516014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM); 1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *); 15186c699258SBarry Smith 1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *); 1520d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *); 15210f5c6efeSJed Brown 1522f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *); 1523f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *); 1524aad739acSMatthew G. Knepley 15252e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec)); 15262e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec)); 15272e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec); 15282e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec)); 1529aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec)); 1530aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec); 1531f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool); 1532d60b7d5cSBarry Smith 1533d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure); 1534