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 15435da068SBarry Smith TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 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 2487497f52SBarry Smith TSType - String with the name of a PETSc `TS` method. 25435da068SBarry Smith 26435da068SBarry Smith Level: beginner 27435da068SBarry Smith 281cc06b55SBarry Smith .seealso: [](integrator_table), [](ch_ts), `TSSetType()`, `TS`, `TSRegister()` 2976bdecfbSBarry Smith J*/ 3019fd82e9SBarry Smith typedef const char *TSType; 319596e0b4SJed Brown #define TSEULER "euler" 329596e0b4SJed Brown #define TSBEULER "beuler" 33e49d4f37SHong Zhang #define TSBASICSYMPLECTIC "basicsymplectic" 349596e0b4SJed Brown #define TSPSEUDO "pseudo" 354d91e141SJed Brown #define TSCN "cn" 369596e0b4SJed Brown #define TSSUNDIALS "sundials" 374d91e141SJed Brown #define TSRK "rk" 389596e0b4SJed Brown #define TSPYTHON "python" 399596e0b4SJed Brown #define TSTHETA "theta" 4088df8a41SLisandro Dalcin #define TSALPHA "alpha" 41818efac9SLisandro Dalcin #define TSALPHA2 "alpha2" 4226d28e4eSEmil Constantinescu #define TSGLLE "glle" 43b6a60446SDebojyoti Ghosh #define TSGLEE "glee" 44ef7bb5aaSJed Brown #define TSSSP "ssp" 458a381b04SJed Brown #define TSARKIMEX "arkimex" 46e27a552bSJed Brown #define TSROSW "rosw" 477d5697caSHong Zhang #define TSEIMEX "eimex" 48abadfa0fSMatthew G. Knepley #define TSMIMEX "mimex" 49211a84d6SLisandro Dalcin #define TSBDF "bdf" 50d249a78fSBarry Smith #define TSRADAU5 "radau5" 514b84eec9SHong Zhang #define TSMPRK "mprk" 5240be0ff1SMatthew G. Knepley #define TSDISCGRAD "discgrad" 53d2567f34SHong Zhang #define TSIRK "irk" 54d27334e2SStefano Zampini #define TSDIRK "dirk" 5540be0ff1SMatthew G. Knepley 56435da068SBarry Smith /*E 5787497f52SBarry Smith TSProblemType - Determines the type of problem this `TS` object is to be used to solve 58435da068SBarry Smith 59195e9b02SBarry Smith Values: 60195e9b02SBarry Smith + `TS_LINEAR` - a linear ODE or DAE 61195e9b02SBarry Smith - `TS_NONLINEAR` - a nonlinear ODE or DAE 62195e9b02SBarry Smith 63435da068SBarry Smith Level: beginner 64435da068SBarry Smith 651cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()` 66435da068SBarry Smith E*/ 679371c9d4SSatish Balay typedef enum { 689371c9d4SSatish Balay TS_LINEAR, 699371c9d4SSatish Balay TS_NONLINEAR 709371c9d4SSatish Balay } TSProblemType; 71818ad0c1SBarry Smith 72b93fea0eSJed Brown /*E 7387497f52SBarry Smith TSEquationType - type of `TS` problem that is solved 74e817cc15SEmil Constantinescu 75195e9b02SBarry Smith Values: 76195e9b02SBarry Smith + `TS_EQ_UNSPECIFIED` - (default) 7716a05f60SBarry 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 7816a05f60SBarry Smith - `TS_EQ_IMPLICIT` - {ODE and DAE index 1, 2, 3, HI} F(t,U,U_t) = 0 79e817cc15SEmil Constantinescu 8095bd0b28SBarry Smith Level: beginner 8195bd0b28SBarry Smith 821cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetEquationType()`, `TSSetEquationType()` 83e817cc15SEmil Constantinescu E*/ 84e817cc15SEmil Constantinescu typedef enum { 85e817cc15SEmil Constantinescu TS_EQ_UNSPECIFIED = -1, 86e817cc15SEmil Constantinescu TS_EQ_EXPLICIT = 0, 87e817cc15SEmil Constantinescu TS_EQ_ODE_EXPLICIT = 1, 88e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 89e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 90e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 91e817cc15SEmil Constantinescu TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 92e817cc15SEmil Constantinescu TS_EQ_IMPLICIT = 1000, 93e817cc15SEmil Constantinescu TS_EQ_ODE_IMPLICIT = 1001, 94e817cc15SEmil Constantinescu TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 95e817cc15SEmil Constantinescu TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 96e817cc15SEmil Constantinescu TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 9719436ca2SJed Brown TS_EQ_DAE_IMPLICIT_INDEXHI = 1500 98e817cc15SEmil Constantinescu } TSEquationType; 99e817cc15SEmil Constantinescu PETSC_EXTERN const char *const *TSEquationTypes; 100e817cc15SEmil Constantinescu 101e817cc15SEmil Constantinescu /*E 10216a05f60SBarry Smith TSConvergedReason - reason a `TS` method has converged (integrated to the requested time) or not 103b93fea0eSJed Brown 104195e9b02SBarry Smith Values: 105195e9b02SBarry Smith + `TS_CONVERGED_ITERATING` - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()` 106195e9b02SBarry Smith . `TS_CONVERGED_TIME` - the final time was reached 107195e9b02SBarry Smith . `TS_CONVERGED_ITS` - the maximum number of iterations (time-steps) was reached prior to the final time 108195e9b02SBarry Smith . `TS_CONVERGED_USER` - user requested termination 109195e9b02SBarry Smith . `TS_CONVERGED_EVENT` - user requested termination on event detection 110195e9b02SBarry Smith . `TS_CONVERGED_PSEUDO_FATOL` - stops when function norm decreased by a set amount, used only for `TSPSEUDO` 111195e9b02SBarry Smith . `TS_CONVERGED_PSEUDO_FRTOL` - stops when function norm decreases below a set amount, used only for `TSPSEUDO` 112195e9b02SBarry Smith . `TS_DIVERGED_NONLINEAR_SOLVE` - too many nonlinear solve failures have occurred 113195e9b02SBarry Smith . `TS_DIVERGED_STEP_REJECTED` - too many steps were rejected 114195e9b02SBarry Smith . `TSFORWARD_DIVERGED_LINEAR_SOLVE` - tangent linear solve failed 115195e9b02SBarry Smith - `TSADJOINT_DIVERGED_LINEAR_SOLVE` - transposed linear solve failed 116195e9b02SBarry Smith 117b93fea0eSJed Brown Level: beginner 118b93fea0eSJed Brown 1191cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()` 120b93fea0eSJed Brown E*/ 121193ac0bcSJed Brown typedef enum { 122193ac0bcSJed Brown TS_CONVERGED_ITERATING = 0, 123193ac0bcSJed Brown TS_CONVERGED_TIME = 1, 124193ac0bcSJed Brown TS_CONVERGED_ITS = 2, 125d6ad946cSShri Abhyankar TS_CONVERGED_USER = 3, 1262b7db910SShri Abhyankar TS_CONVERGED_EVENT = 4, 1273118ae5eSBarry Smith TS_CONVERGED_PSEUDO_FATOL = 5, 1283118ae5eSBarry Smith TS_CONVERGED_PSEUDO_FRTOL = 6, 129193ac0bcSJed Brown TS_DIVERGED_NONLINEAR_SOLVE = -1, 130b5b4017aSHong Zhang TS_DIVERGED_STEP_REJECTED = -2, 13105c9950eSHong Zhang TSFORWARD_DIVERGED_LINEAR_SOLVE = -3, 13205c9950eSHong Zhang TSADJOINT_DIVERGED_LINEAR_SOLVE = -4 133193ac0bcSJed Brown } TSConvergedReason; 134014dd563SJed Brown PETSC_EXTERN const char *const *TSConvergedReasons; 1351957e957SBarry Smith 136b93fea0eSJed Brown /*MC 13787497f52SBarry Smith TS_CONVERGED_ITERATING - this only occurs if `TSGetConvergedReason()` is called during the `TSSolve()` 138b93fea0eSJed Brown 139b93fea0eSJed Brown Level: beginner 140b93fea0eSJed Brown 1411cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()` 142b93fea0eSJed Brown M*/ 143b93fea0eSJed Brown 144b93fea0eSJed Brown /*MC 145b93fea0eSJed Brown TS_CONVERGED_TIME - the final time was reached 146b93fea0eSJed Brown 147b93fea0eSJed Brown Level: beginner 148b93fea0eSJed Brown 1491cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxTime()`, `TSGetMaxTime()`, `TSGetSolveTime()` 150b93fea0eSJed Brown M*/ 151b93fea0eSJed Brown 152b93fea0eSJed Brown /*MC 1531957e957SBarry Smith TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time 154b93fea0eSJed Brown 155b93fea0eSJed Brown Level: beginner 156b93fea0eSJed Brown 1571cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxSteps()`, `TSGetMaxSteps()` 158b93fea0eSJed Brown M*/ 1591957e957SBarry Smith 160d6ad946cSShri Abhyankar /*MC 161d6ad946cSShri Abhyankar TS_CONVERGED_USER - user requested termination 162d6ad946cSShri Abhyankar 163d6ad946cSShri Abhyankar Level: beginner 164d6ad946cSShri Abhyankar 1651cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()` 166b93fea0eSJed Brown M*/ 167b93fea0eSJed Brown 168b93fea0eSJed Brown /*MC 169d76fc68cSShri Abhyankar TS_CONVERGED_EVENT - user requested termination on event detection 170d76fc68cSShri Abhyankar 171d76fc68cSShri Abhyankar Level: beginner 172d76fc68cSShri Abhyankar 1731cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()` 174d76fc68cSShri Abhyankar M*/ 175d76fc68cSShri Abhyankar 176d76fc68cSShri Abhyankar /*MC 17787497f52SBarry Smith TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for `TSPSEUDO` 1783118ae5eSBarry Smith 1793c7db156SBarry Smith Options Database Key: 180d5b43468SJose E. Roman . -ts_pseudo_frtol <rtol> - use specified rtol 1813118ae5eSBarry Smith 182195e9b02SBarry Smith Level: beginner 183195e9b02SBarry Smith 1841cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FATOL` 1853118ae5eSBarry Smith M*/ 1863118ae5eSBarry Smith 1873118ae5eSBarry Smith /*MC 18887497f52SBarry Smith TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for `TSPSEUDO` 1893118ae5eSBarry Smith 1903c7db156SBarry Smith Options Database Key: 19167b8a455SSatish Balay . -ts_pseudo_fatol <atol> - use specified atol 1923118ae5eSBarry Smith 193195e9b02SBarry Smith Level: beginner 194195e9b02SBarry Smith 1951cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSSetConvergedReason()`, `TS_CONVERGED_PSEUDO_FRTOL` 1963118ae5eSBarry Smith M*/ 1973118ae5eSBarry Smith 1983118ae5eSBarry Smith /*MC 199b93fea0eSJed Brown TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 200b93fea0eSJed Brown 201b93fea0eSJed Brown Level: beginner 202b93fea0eSJed Brown 203195e9b02SBarry Smith Note: 204195e9b02SBarry Smith See `TSSetMaxSNESFailures()` for how to allow more nonlinear solver failures. 2051957e957SBarry Smith 2061cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSGetSNES()`, `SNESGetConvergedReason()`, `TSSetMaxSNESFailures()` 207b93fea0eSJed Brown M*/ 208b93fea0eSJed Brown 209b93fea0eSJed Brown /*MC 210b93fea0eSJed Brown TS_DIVERGED_STEP_REJECTED - too many steps were rejected 211b93fea0eSJed Brown 212b93fea0eSJed Brown Level: beginner 213b93fea0eSJed Brown 21495452b02SPatrick Sanan Notes: 215195e9b02SBarry Smith See `TSSetMaxStepRejections()` for how to allow more step rejections. 2161957e957SBarry Smith 2171cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSGetConvergedReason()`, `TSGetAdapt()`, `TSSetMaxStepRejections()` 218b93fea0eSJed Brown M*/ 219b93fea0eSJed Brown 220d6ad946cSShri Abhyankar /*E 221d6ad946cSShri Abhyankar TSExactFinalTimeOption - option for handling of final time step 222d6ad946cSShri Abhyankar 223195e9b02SBarry Smith Values: 22416a05f60SBarry Smith + `TS_EXACTFINALTIME_STEPOVER` - Don't do anything if requested final time is exceeded 225195e9b02SBarry Smith . `TS_EXACTFINALTIME_INTERPOLATE` - Interpolate back to final time 22616a05f60SBarry Smith - `TS_EXACTFINALTIME_MATCHSTEP` - Adapt final time step to match the final time requested 227195e9b02SBarry Smith 228d6ad946cSShri Abhyankar Level: beginner 229d6ad946cSShri Abhyankar 2301cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetConvergedReason()`, `TSSetExactFinalTime()`, `TSGetExactFinalTime()` 231d6ad946cSShri Abhyankar E*/ 2329371c9d4SSatish Balay typedef enum { 2339371c9d4SSatish Balay TS_EXACTFINALTIME_UNSPECIFIED = 0, 2349371c9d4SSatish Balay TS_EXACTFINALTIME_STEPOVER = 1, 2359371c9d4SSatish Balay TS_EXACTFINALTIME_INTERPOLATE = 2, 2369371c9d4SSatish Balay TS_EXACTFINALTIME_MATCHSTEP = 3 2379371c9d4SSatish Balay } TSExactFinalTimeOption; 238d6ad946cSShri Abhyankar PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 239d6ad946cSShri Abhyankar 240000e7ae3SMatthew Knepley /* Logging support */ 241014dd563SJed Brown PETSC_EXTERN PetscClassId TS_CLASSID; 242d74926cbSBarry Smith PETSC_EXTERN PetscClassId DMTS_CLASSID; 2431b9b13dfSLisandro Dalcin PETSC_EXTERN PetscClassId TSADAPT_CLASSID; 244000e7ae3SMatthew Knepley 245607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 246560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSFinalizePackage(void); 2478ba1e511SMatthew Knepley 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm, TS *); 249baa10174SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSClone(TS, TS *); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSDestroy(TS *); 251818ad0c1SBarry Smith 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetProblemType(TS, TSProblemType); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetProblemType(TS, TSProblemType *); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitor(TS, PetscInt, PetscReal, Vec); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, void *), void *, PetscErrorCode (*)(void **)); 256721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *)); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 258818ad0c1SBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS, const char[]); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS, const char[]); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS, const char *[]); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetUp(TS); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSReset(TS); 265818ad0c1SBarry Smith 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetSolution(TS, Vec); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolution(TS, Vec *); 268818ad0c1SBarry Smith 269efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2SetSolution(TS, Vec, Vec); 270efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TS2GetSolution(TS, Vec *, Vec *); 27157df6a1bSDebojyoti Ghosh 272642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetSolutionComponents(TS, PetscInt *, Vec *); 273642f3702SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetAuxSolution(TS, Vec *); 2740a01e1b2SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetTimeError(TS, PetscInt, Vec *); 27557df6a1bSDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetTimeError(TS, Vec); 2766e2e232bSDebojyoti Ghosh 277a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetRHSJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *); 278cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetRHSJacobianP(TS, Mat *, PetscErrorCode (**)(TS, PetscReal, Vec, Mat, void *), void **); 279a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianP(TS, PetscReal, Vec, Mat); 28033f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSSetIJacobianP(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *); 28133f72c5dSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIJacobianP(TS, PetscReal, Vec, Vec, PetscReal, Mat, PetscBool); 282edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobianP()", ) PetscErrorCode TSComputeDRDPFunction(TS, PetscReal, Vec, Vec *); 283edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS() then TSComputeRHSJacobian()", ) PetscErrorCode TSComputeDRDUFunction(TS, PetscReal, Vec, Vec *); 2849371c9d4SSatish 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 *); 285b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 286b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 287b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 288b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeIHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 2899371c9d4SSatish 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 *); 290b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 291b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 292b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS, PetscReal, Vec, Vec *, Vec, Vec *); 293b1e111ebSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS, PetscReal, Vec, Vec *, Vec, Vec *); 294b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostHessianProducts(TS, PetscInt, Vec *, Vec *, Vec); 295b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS, PetscInt *, Vec **, Vec **, Vec *); 296ba0675f6SHong Zhang PETSC_EXTERN PetscErrorCode TSComputeSNESJacobian(TS, Vec, Mat, Mat); 297a05bf03eSHong Zhang 298bc952696SBarry Smith /*S 299e91eccc2SStefano Zampini TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step) 300bc952696SBarry Smith 301bc952696SBarry Smith Level: advanced 302bc952696SBarry Smith 3031cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectorySetType()`, `TSTrajectoryDestroy()`, `TSTrajectoryReset()` 304bc952696SBarry Smith S*/ 305bc952696SBarry Smith typedef struct _p_TSTrajectory *TSTrajectory; 306bc952696SBarry Smith 307bc952696SBarry Smith /*J 308195e9b02SBarry Smith TSTrajectoryType - String with the name of a PETSc `TS` trajectory storage method 309bc952696SBarry Smith 310bc952696SBarry Smith Level: intermediate 311bc952696SBarry Smith 3121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetSaveTrajectory()`, `TSTrajectoryCreate()`, `TSTrajectoryDestroy()` 313bc952696SBarry Smith J*/ 314bc952696SBarry Smith typedef const char *TSTrajectoryType; 315bc952696SBarry Smith #define TSTRAJECTORYBASIC "basic" 3161c8c567eSBarry Smith #define TSTRAJECTORYSINGLEFILE "singlefile" 3179a53571cSHong Zhang #define TSTRAJECTORYMEMORY "memory" 3182b043167SHong Zhang #define TSTRAJECTORYVISUALIZATION "visualization" 319bc952696SBarry Smith 320bc952696SBarry Smith PETSC_EXTERN PetscFunctionList TSTrajectoryList; 321bc952696SBarry Smith PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID; 322bc952696SBarry Smith PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled; 323bc952696SBarry Smith 324bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS); 3252d29f1f2SStefano Zampini PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS); 32667a3cfb0SHong Zhang PETSC_EXTERN PetscErrorCode TSRemoveTrajectory(TS); 327bc952696SBarry Smith 328bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm, TSTrajectory *); 3299a992471SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory); 330bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory *); 331560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory, PetscViewer); 332fd9d3c67SJed Brown PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory, TS, TSTrajectoryType); 333881c1a9bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGetType(TSTrajectory, TS, TSTrajectoryType *); 334bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory, TS, PetscInt, PetscReal, Vec); 335c679fc15SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory, TS, PetscInt, PetscReal *); 336fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory, TS, PetscInt, PetscReal *, Vec, Vec); 337fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory, TS, PetscReal, Vec *, Vec *); 338fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory, PetscInt *); 339fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory, Vec *, Vec *); 340972caf09SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory, TS); 341560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[], PetscErrorCode (*)(TSTrajectory, TS)); 342bc952696SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void); 34368bece0bSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory, TS); 3449ffb3502SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetUseHistory(TSTrajectory, PetscBool); 3452bee684fSHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory, PetscBool); 34678fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory, const char *const *); 34708347785SBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 348fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory, PetscBool); 349fe8322adSStefano Zampini PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory, PetscBool *); 35064fc91eeSBarry Smith PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory, PetscBool); 35164e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory, const char[]); 35264e38db7SHong Zhang PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory, const char[]); 35378fbdcc8SBarry Smith PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS, TSTrajectory *); 354bc952696SBarry Smith 3554bf303faSJacob Faibussowitsch typedef enum { 3564bf303faSJacob Faibussowitsch TJ_REVOLVE, 3574bf303faSJacob Faibussowitsch TJ_CAMS, 3584bf303faSJacob Faibussowitsch TJ_PETSC 3594bf303faSJacob Faibussowitsch } TSTrajectoryMemoryType; 3604bf303faSJacob Faibussowitsch PETSC_EXTERN const char *const TSTrajectoryMemoryTypes[]; 3614bf303faSJacob Faibussowitsch 3624bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectoryMemorySetType(TSTrajectory, TSTrajectoryMemoryType); 3634bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsRAM(TSTrajectory, PetscInt); 3644bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxCpsDisk(TSTrajectory, PetscInt); 3654bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsRAM(TSTrajectory, PetscInt); 3664bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSTrajectorySetMaxUnitsDisk(TSTrajectory, PetscInt); 3674bf303faSJacob Faibussowitsch 368dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS, PetscInt, Vec *, Vec *); 369dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS, PetscInt *, Vec **, Vec **); 370edd03b47SJacob 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 *); 371dfb21088SHong Zhang PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS, Vec *); 372022c081aSHong Zhang PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS, PetscReal, Vec, Vec); 373cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSCreateQuadratureTS(TS, PetscBool, TS *); 374cd4cee2dSHong Zhang PETSC_EXTERN PetscErrorCode TSGetQuadratureTS(TS, PetscBool *, TS *); 375d4aa0a58SBarry Smith 376dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(TS, PetscOptionItems *); 377a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *); 378a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *, PetscErrorCode (*)(void **)); 379a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS); 380a05bf03eSHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS, const char[], const char[], const char[], PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*)(TS, PetscViewerAndFormat *)); 381a05bf03eSHong Zhang 382edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSSetRHSJacobianP()", ) PetscErrorCode TSAdjointSetRHSJacobian(TS, Mat, PetscErrorCode (*)(TS, PetscReal, Vec, Mat, void *), void *); 383edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSComputeRHSJacobianP()", ) PetscErrorCode TSAdjointComputeRHSJacobian(TS, PetscReal, Vec, Mat); 384edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDPFunction(TS, PetscReal, Vec, Vec *); 385edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 5, 0, "TSGetQuadratureTS()", ) PetscErrorCode TSAdjointComputeDRDYFunction(TS, PetscReal, Vec, Vec *); 3862c39e106SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS); 38773b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS, PetscInt); 388d4aa0a58SBarry Smith 38973b18844SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 390f6a906c0SBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 391ece46509SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointReset(TS); 392999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS); 393ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointSetForward(TS, Mat); 394ecf68647SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointResetForward(TS); 395715f1b00SHong Zhang 39613b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS, PetscInt, Mat); 39713b1b0ffSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS, PetscInt *, Mat *); 398edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSCreateQuadratureTS()", ) PetscErrorCode TSForwardSetIntegralGradients(TS, PetscInt, Vec *); 399edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 12, 0, "TSForwardGetSensitivities()", ) PetscErrorCode TSForwardGetIntegralGradients(TS, PetscInt *, Vec **); 400715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS); 4017adebddeSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardReset(TS); 402999a3721SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS); 403715f1b00SHong Zhang PETSC_EXTERN PetscErrorCode TSForwardStep(TS); 404b5b4017aSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS, Mat); 405cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS, PetscInt *, Mat *[]); 406c235aa8dSHong Zhang 407618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS, PetscInt); 408618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS, PetscInt *); 409618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS, PetscReal); 410618ce8baSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS, PetscReal *); 41149354f04SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS, TSExactFinalTimeOption); 412f6953c82SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS, TSExactFinalTimeOption *); 4134a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSSetTimeSpan(TS, PetscInt, PetscReal *); 4144a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpan(TS, PetscInt *, const PetscReal **); 4154a658b32SHong Zhang PETSC_EXTERN PetscErrorCode TSGetTimeSpanSolutions(TS, PetscInt *, Vec **); 416818ad0c1SBarry Smith 417edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetTime()", ) PetscErrorCode TSSetInitialTimeStep(TS, PetscReal, PetscReal); 418edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSSetMax()", ) PetscErrorCode TSSetDuration(TS, PetscInt, PetscReal); 419edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetMax()", ) PetscErrorCode TSGetDuration(TS, PetscInt *, PetscReal *); 420edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTimeStepNumber(TS, PetscInt *); 421edd03b47SJacob Faibussowitsch PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 8, 0, "TSGetStepNumber()", ) PetscErrorCode TSGetTotalSteps(TS, PetscInt *); 42219eac22cSLisandro Dalcin 423721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 424cc9c3a59SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 42583a4ac43SBarry Smith 42683a4ac43SBarry Smith typedef struct _n_TSMonitorDrawCtx *TSMonitorDrawCtx; 42783a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorDrawCtx *); 42883a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx *); 42983a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS, PetscInt, PetscReal, Vec, void *); 4302d5ee99bSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS, PetscInt, PetscReal, Vec, void *); 43183a4ac43SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS, PetscInt, PetscReal, Vec, void *); 4320ed3bfb6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS, PetscInt, PetscReal, Vec, void *); 43383a4ac43SBarry Smith 434721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *); 4359110b2e7SHong Zhang PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *); 4369110b2e7SHong Zhang 437721cd6eeSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 438014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS, PetscInt, PetscReal, Vec, void *); 439014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void *); 440fb1732b5SBarry Smith 441014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSStep(TS); 4427cbde773SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS, NormType, PetscInt *, PetscReal *); 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS, PetscInt, Vec, PetscBool *); 444cc708dedSBarry Smith PETSC_EXTERN PetscErrorCode TSSolve(TS, Vec); 445e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGetEquationType(TS, TSEquationType *); 446e817cc15SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSSetEquationType(TS, TSEquationType); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS, TSConvergedReason *); 448d6ad946cSShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS, TSConvergedReason); 449487e0bb9SJed Brown PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS, PetscReal *); 4505ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS, PetscInt *); 4515ef26d82SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS, PetscInt *); 452cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS, PetscInt *); 453cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS, PetscInt); 454cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS, PetscInt *); 455cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS, PetscInt); 456cef5090cSJed Brown PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS, PetscBool); 457dcb233daSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRestartStep(TS); 45824655328SShri PETSC_EXTERN PetscErrorCode TSRollBack(TS); 459d2daff3dSHong Zhang 460cc4f23bcSHong Zhang PETSC_EXTERN PetscErrorCode TSGetStages(TS, PetscInt *, Vec *[]); 461818ad0c1SBarry Smith 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetTime(TS, PetscReal *); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTime(TS, PetscReal); 464e5e524a1SHong Zhang PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS, PetscReal *); 46580275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS, PetscReal *); 46680275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS, PetscReal); 46780275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS, PetscInt *); 46880275a0aSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS, PetscInt); 469818ad0c1SBarry Smith 47014d0ab18SJacob Faibussowitsch /*S 47114d0ab18SJacob Faibussowitsch TSRHSFunction - An `TS` right-hand-side evaluation function 47214d0ab18SJacob Faibussowitsch 47314d0ab18SJacob Faibussowitsch Calling Sequence: 47414d0ab18SJacob Faibussowitsch + ts - timestep context 47514d0ab18SJacob Faibussowitsch . t - current time 47614d0ab18SJacob Faibussowitsch . u - input vector 47714d0ab18SJacob Faibussowitsch . F - function vector 47814d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 47914d0ab18SJacob Faibussowitsch 48014d0ab18SJacob Faibussowitsch Level: beginner 48114d0ab18SJacob Faibussowitsch 48214d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `TSIFunction()`, 48314d0ab18SJacob Faibussowitsch `TSIJacobian()`, `TSRHSJacobian()` 48414d0ab18SJacob Faibussowitsch S*/ 48514d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS ts, PetscReal t, Vec u, Vec F, void *ctx); 48614d0ab18SJacob Faibussowitsch 48714d0ab18SJacob Faibussowitsch /*S 48814d0ab18SJacob Faibussowitsch TSRHSJacobian - A `TS` right-hand-side Jacobian evaluation function 48914d0ab18SJacob Faibussowitsch 49014d0ab18SJacob Faibussowitsch Calling Sequence: 49114d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 49214d0ab18SJacob Faibussowitsch . t - current time 49314d0ab18SJacob Faibussowitsch . u - input vector 49414d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian matrix 49514d0ab18SJacob Faibussowitsch . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat) 49614d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for matrix evaluation routine 49714d0ab18SJacob Faibussowitsch 49814d0ab18SJacob Faibussowitsch Level: beginner 49914d0ab18SJacob Faibussowitsch 50014d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobian()`, `DMTSSetRHSJacobian()`, `TSRHSFunction()`, 50114d0ab18SJacob Faibussowitsch `TSIFunction()`, `TSIJacobian()` 50214d0ab18SJacob Faibussowitsch S*/ 50314d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS ts, PetscReal t, Vec u, Mat Amat, Mat Pmat, void *ctx); 50414d0ab18SJacob Faibussowitsch 50514d0ab18SJacob Faibussowitsch /*S 50614d0ab18SJacob Faibussowitsch TSRHSJacobianP - A function that computes the Jacobian of G w.r.t. the parameters P where U_t 50714d0ab18SJacob Faibussowitsch = G(U,P,t), as well as the location to store the matrix. 50814d0ab18SJacob Faibussowitsch 50914d0ab18SJacob Faibussowitsch Calling Sequence: 51014d0ab18SJacob Faibussowitsch + ts - the `TS` context 51114d0ab18SJacob Faibussowitsch . t - current timestep 51214d0ab18SJacob Faibussowitsch . U - input vector (current ODE solution) 51314d0ab18SJacob Faibussowitsch . A - output matrix 51414d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 51514d0ab18SJacob Faibussowitsch 51614d0ab18SJacob Faibussowitsch Level: beginner 51714d0ab18SJacob Faibussowitsch 51814d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetRHSJacobianP()` 51914d0ab18SJacob Faibussowitsch S*/ 52014d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobianP)(TS ts, PetscReal t, Vec U, Mat A, void *ctx); 52114d0ab18SJacob Faibussowitsch 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS, Vec, TSRHSFunction, void *); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS, Vec *, TSRHSFunction *, void **); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS, Mat, Mat, TSRHSJacobian, void *); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS, Mat *, Mat *, TSRHSJacobian *, void **); 526e1244c69SJed Brown PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS, PetscBool); 527818ad0c1SBarry Smith 52814d0ab18SJacob Faibussowitsch /*S 52914d0ab18SJacob Faibussowitsch TSSolutionFunction - A `TS` solution evaluation function 53014d0ab18SJacob Faibussowitsch 53114d0ab18SJacob Faibussowitsch Calling Sequence: 53214d0ab18SJacob Faibussowitsch + ts - timestep context 53314d0ab18SJacob Faibussowitsch . t - current time 53414d0ab18SJacob Faibussowitsch . u - output vector 53514d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 53614d0ab18SJacob Faibussowitsch 53714d0ab18SJacob Faibussowitsch Level: advanced 53814d0ab18SJacob Faibussowitsch 53914d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetSolutionFunction()`, `DMTSSetSolutionFunction()` 54014d0ab18SJacob Faibussowitsch S*/ 54114d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS ts, PetscReal t, Vec u, void *ctx); 54214d0ab18SJacob Faibussowitsch 543ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS, TSSolutionFunction, void *); 54414d0ab18SJacob Faibussowitsch 54514d0ab18SJacob Faibussowitsch /*S 54614d0ab18SJacob Faibussowitsch TSForcingFunction - A `TS` forcing function evaluation function 54714d0ab18SJacob Faibussowitsch 54814d0ab18SJacob Faibussowitsch Calling Sequence: 54914d0ab18SJacob Faibussowitsch + ts - timestep context 55014d0ab18SJacob Faibussowitsch . t - current time 55114d0ab18SJacob Faibussowitsch . f - output vector 55214d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 55314d0ab18SJacob Faibussowitsch 55414d0ab18SJacob Faibussowitsch Level: advanced 55514d0ab18SJacob Faibussowitsch 55614d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetForcingFunction()`, `DMTSSetForcingFunction()` 55714d0ab18SJacob Faibussowitsch S*/ 55814d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS ts, PetscReal t, Vec f, void *ctx); 55914d0ab18SJacob Faibussowitsch 560d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS, TSForcingFunction, void *); 561ef20d060SBarry Smith 56214d0ab18SJacob Faibussowitsch /*S 56314d0ab18SJacob Faibussowitsch TSIFunction - A `TS` implicit function evaluation function 56414d0ab18SJacob Faibussowitsch 56514d0ab18SJacob Faibussowitsch Calling Sequence: 56614d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 56714d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 56814d0ab18SJacob Faibussowitsch . U - state vector 56914d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 57014d0ab18SJacob Faibussowitsch . F - function vector 57114d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for function 57214d0ab18SJacob Faibussowitsch 57314d0ab18SJacob Faibussowitsch Level: beginner 57414d0ab18SJacob Faibussowitsch 57514d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetIFunction()`, `DMTSSetIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()` 57614d0ab18SJacob Faibussowitsch S*/ 57714d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS ts, PetscReal t, Vec U, Vec U_tt, Vec F, void *ctx); 57814d0ab18SJacob Faibussowitsch 57914d0ab18SJacob Faibussowitsch /*S 58014d0ab18SJacob Faibussowitsch TSIJacobian - A `TS` Jacobian evaluation function 58114d0ab18SJacob Faibussowitsch 58214d0ab18SJacob Faibussowitsch Calling Sequence: 58314d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 58414d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 58514d0ab18SJacob Faibussowitsch . U - state vector 58614d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 58714d0ab18SJacob Faibussowitsch . a - shift 58814d0ab18SJacob Faibussowitsch . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 58914d0ab18SJacob Faibussowitsch . Pmat - matrix used for constructing preconditioner, usually the same as `Amat` 59014d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for Jacobian evaluation routine 59114d0ab18SJacob Faibussowitsch 59214d0ab18SJacob Faibussowitsch Level: beginner 59314d0ab18SJacob Faibussowitsch 59414d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TSSetIJacobian()`, `DMTSSetIJacobian()`, `TSIFunction()`, `TSRHSFunction()`, `TSRHSJacobian()` 59514d0ab18SJacob Faibussowitsch S*/ 59614d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS ts, PetscReal t, Vec U, Vec U_t, PetscReal a, Mat Amat, Mat Pmat, void *ctx); 59714d0ab18SJacob Faibussowitsch 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetIFunction(TS, Vec, TSIFunction, void *); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetIFunction(TS, Vec *, TSIFunction *, void **); 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS, Mat, Mat, TSIJacobian, void *); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS, Mat *, Mat *, TSIJacobian *, void **); 602316643e7SJed Brown 60314d0ab18SJacob Faibussowitsch /*S 60414d0ab18SJacob Faibussowitsch TSI2Function - A `TS` implicit function evaluation function for 2nd order systems 60514d0ab18SJacob Faibussowitsch 60614d0ab18SJacob Faibussowitsch Calling Sequence: 60714d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 60814d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 60914d0ab18SJacob Faibussowitsch . U - state vector 61014d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 61114d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector 61214d0ab18SJacob Faibussowitsch . F - function vector 61314d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for matrix evaluation routine (may be `NULL`) 61414d0ab18SJacob Faibussowitsch 61514d0ab18SJacob Faibussowitsch Level: advanced 61614d0ab18SJacob Faibussowitsch 61714d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetI2Function()`, `DMTSSetI2Function()` 61814d0ab18SJacob Faibussowitsch S*/ 61914d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, Vec F, void *ctx); 62014d0ab18SJacob Faibussowitsch 62114d0ab18SJacob Faibussowitsch /*S 62214d0ab18SJacob Faibussowitsch TSI2Jacobian - A `TS` implicit Jacobian evaluation function for 2nd order systems 62314d0ab18SJacob Faibussowitsch 62414d0ab18SJacob Faibussowitsch Calling Sequence: 62514d0ab18SJacob Faibussowitsch + ts - the `TS` context obtained from `TSCreate()` 62614d0ab18SJacob Faibussowitsch . t - time at step/stage being solved 62714d0ab18SJacob Faibussowitsch . U - state vector 62814d0ab18SJacob Faibussowitsch . U_t - time derivative of state vector 62914d0ab18SJacob Faibussowitsch . U_tt - second time derivative of state vector 63014d0ab18SJacob Faibussowitsch . v - shift for U_t 63114d0ab18SJacob Faibussowitsch . a - shift for U_tt 63214d0ab18SJacob 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 63314d0ab18SJacob Faibussowitsch . jac - matrix from which to construct the preconditioner, may be same as `J` 63414d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined context for matrix evaluation routine 63514d0ab18SJacob Faibussowitsch 63614d0ab18SJacob Faibussowitsch Level: advanced 63714d0ab18SJacob Faibussowitsch 63814d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetI2Jacobian()`, `DMTSSetI2Jacobian()`, `TSIFunction()`, `TSIJacobian()`, `TSRHSFunction()`, `TSRHSJacobian()` 63914d0ab18SJacob Faibussowitsch S*/ 64014d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS ts, PetscReal t, Vec U, Vec U_t, Vec U_tt, PetscReal v, PetscReal a, Mat J, Mat Jac, void *ctx); 64114d0ab18SJacob Faibussowitsch 642efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetI2Function(TS, Vec, TSI2Function, void *); 643efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetI2Function(TS, Vec *, TSI2Function *, void **); 644efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS, Mat, Mat, TSI2Jacobian, void *); 645efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS, Mat *, Mat *, TSI2Jacobian *, void **); 646efe9872eSLisandro Dalcin 647545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS, const char[], IS); 648545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS, const char[], IS *); 649545aaa6fSHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS, const char[], Vec, TSRHSFunction, void *); 6501d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS, const char[], TS *); 6511d06f6b3SHong Zhang PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS, PetscInt *, TS *[]); 6520fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool); 6530fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool *); 654d9194312SHong Zhang 655014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS, PetscReal, Vec, Vec, void *); 656d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS, PetscReal, Vec, Mat, Mat, void *); 657014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS, PetscReal, Vec, Vec, Vec, void *); 658d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 6591c8a10a0SJed Brown PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS, PetscReal, Vec); 6609b7cd975SBarry Smith PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS, PetscReal, Vec); 661847ff0e1SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 662d7cfae9bSHong Zhang PETSC_EXTERN PetscErrorCode TSPruneIJacobianColor(TS, Mat, Mat); 663e34be4c2SBarry Smith 664014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 665b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS, PetscReal)); 6669be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS, PetscReal, PetscInt, Vec *)); 667c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode (*)(TS)); 668014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 6696bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSSetResize(TS, PetscErrorCode (*)(TS, PetscInt, PetscReal, Vec, PetscBool *, void *), PetscErrorCode (*)(TS, PetscInt, Vec[], Vec[], void *), void *); 670014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPreStep(TS); 671b8123daeSJed Brown PETSC_EXTERN PetscErrorCode TSPreStage(TS, PetscReal); 6729be3e283SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSPostStage(TS, PetscReal, PetscInt, Vec *); 673c688d042SShri Abhyankar PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS); 674014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPostStep(TS); 6756bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResize(TS); 6766bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRetrieveVec(TS, const char *, Vec *); 6776bd3a4fdSStefano Zampini PETSC_EXTERN PetscErrorCode TSResizeRegisterVec(TS, const char *, Vec); 6786bd3a4fdSStefano Zampini 679014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSInterpolate(TS, PetscReal, Vec); 680014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetTolerances(TS, PetscReal, Vec, PetscReal, Vec); 681c5033834SJed Brown PETSC_EXTERN PetscErrorCode TSGetTolerances(TS, PetscReal *, Vec *, PetscReal *, Vec *); 6827453f775SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *); 6838a175baeSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS, Vec, Vec, Vec, NormType, PetscReal *, PetscReal *, PetscReal *); 684014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS, PetscReal); 685014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS, PetscReal *); 686d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS, PetscReal, Vec, PetscBool *)); 687d183316bSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS, PetscReal, Vec, PetscBool *); 688000e7ae3SMatthew Knepley 689014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS, PetscErrorCode (*)(TS, PetscReal *, void *), void *); 6908d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS, PetscReal *, void *); 691014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS, PetscReal *); 692014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS, PetscReal); 693014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS, PetscErrorCode (*)(TS, Vec, void *, PetscReal *, PetscBool *), void *); 6948d359177SBarry Smith PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS, Vec, void *, PetscReal *, PetscBool *); 695014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS, Vec, PetscReal *, PetscBool *); 696014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS, PetscReal); 697014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 69821c89e3eSBarry Smith 699014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSPythonSetType(TS, const char[]); 700ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode TSPythonGetType(TS, const char *[]); 7011d6018f0SLisandro Dalcin 702014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS, PetscReal, Vec, Vec); 703d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS, PetscReal, Vec, Mat, Mat); 704014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS, PetscReal, Vec, Vec, Vec, PetscBool); 705d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, PetscBool); 706efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS, PetscReal, Vec, Vec, Vec, Vec); 707efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS, PetscReal, Vec, Vec, Vec, PetscReal, PetscReal, Mat, Mat); 708f9c1d6abSBarry Smith PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS, PetscReal, PetscReal, PetscReal *, PetscReal *); 709818ad0c1SBarry Smith 710014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS, Vec, Vec); 711d6ebe24aSShri Abhyankar 712967bb25aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 71324989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM, TSRHSFunction, void *); 71424989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM, TSRHSFunction *, void **); 715800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionContextDestroy(DM, PetscErrorCode (*)(void *)); 71624989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM, TSRHSJacobian, void *); 71724989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM, TSRHSJacobian *, void **); 718800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobianContextDestroy(DM, PetscErrorCode (*)(void *)); 71924989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM, TSIFunction, void *); 72024989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM, TSIFunction *, void **); 721800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIFunctionContextDestroy(DM, PetscErrorCode (*)(void *)); 72224989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM, TSIJacobian, void *); 72324989b8cSPeter Brune PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM, TSIJacobian *, void **); 724800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetIJacobianContextDestroy(DM, PetscErrorCode (*)(void *)); 725efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM, TSI2Function, void *); 726efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM, TSI2Function *, void **); 727800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2FunctionContextDestroy(DM, PetscErrorCode (*)(void *)); 728efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM, TSI2Jacobian, void *); 729efe9872eSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM, TSI2Jacobian *, void **); 730800f99ffSJeremy L Thompson PETSC_EXTERN PetscErrorCode DMTSSetI2JacobianContextDestroy(DM, PetscErrorCode (*)(void *)); 731efe9872eSLisandro Dalcin 73214d0ab18SJacob Faibussowitsch /*S 73314d0ab18SJacob Faibussowitsch TSTransientVariable - A function to transform from state to transient variables 73414d0ab18SJacob Faibussowitsch 73514d0ab18SJacob Faibussowitsch Calling Sequence: 73614d0ab18SJacob Faibussowitsch + ts - timestep context 73714d0ab18SJacob Faibussowitsch . p - input vector (primitive form) 73814d0ab18SJacob Faibussowitsch . c - output vector, transient variables (conservative form) 73914d0ab18SJacob Faibussowitsch - ctx - [optional] user-defined function context 74014d0ab18SJacob Faibussowitsch 74114d0ab18SJacob Faibussowitsch Level: advanced 74214d0ab18SJacob Faibussowitsch 74314d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TS`, `TSSetTransientVariable()`, `DMTSSetTransientVariable()` 74414d0ab18SJacob Faibussowitsch S*/ 74514d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSTransientVariable)(TS ts, Vec p, Vec c, void *ctx); 74614d0ab18SJacob Faibussowitsch 747438f35afSJed Brown PETSC_EXTERN PetscErrorCode TSSetTransientVariable(TS, TSTransientVariable, void *); 748e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode DMTSSetTransientVariable(DM, TSTransientVariable, void *); 749e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode DMTSGetTransientVariable(DM, TSTransientVariable *, void *); 750e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSComputeTransientVariable(TS, Vec, Vec); 751e3c11fc1SJed Brown PETSC_EXTERN PetscErrorCode TSHasTransientVariable(TS, PetscBool *); 752e3c11fc1SJed Brown 753ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM, TSSolutionFunction, void *); 754ef20d060SBarry Smith PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM, TSSolutionFunction *, void **); 755d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM, TSForcingFunction, void *); 756d56366bfSLisandro Dalcin PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM, TSForcingFunction *, void **); 757e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckResidual(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscReal *); 758e17bd7bbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckJacobian(TS, DM, PetscReal, Vec, Vec, PetscReal, PetscBool *, PetscReal *); 7597f96f943SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec); 7609b7cd975SBarry Smith 7610c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, Vec, void *), void **); 762d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, Vec, void *), void *); 7630c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetIJacobianLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **); 764d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *); 7650c9409e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSGetRHSFunctionLocal(DM, PetscErrorCode (**)(DM, PetscReal, Vec, Vec, void *), void **); 766d433e6cbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 767b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrix(DM); 768b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM); 769b4937a87SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMTSDestroyRHSMassMatrix(DM); 7706c6b9e74SPeter Brune 7716c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer)); 7726c6b9e74SPeter Brune PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM, PetscErrorCode (*)(void *, PetscViewer), PetscErrorCode (*)(void **, PetscViewer)); 7736c6b9e74SPeter Brune 77414d0ab18SJacob Faibussowitsch /*S 77514d0ab18SJacob Faibussowitsch DMDATSRHSFunctionLocal - A local `TS` right hand side residual evaluation function for use with `DMDA` 7766c6b9e74SPeter Brune 77714d0ab18SJacob Faibussowitsch Calling Sequence: 77814d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 77914d0ab18SJacob Faibussowitsch . t - time at which to evaluate residual 78014d0ab18SJacob Faibussowitsch . x - array of local state information 78114d0ab18SJacob Faibussowitsch . f - output array of local residual information 78214d0ab18SJacob Faibussowitsch - ctx - optional user context 78314d0ab18SJacob Faibussowitsch 78414d0ab18SJacob Faibussowitsch Level: beginner 78514d0ab18SJacob Faibussowitsch 78614d0ab18SJacob Faibussowitsch .seealso: `DMDA`, `DMDATSSetRHSFunctionLocal()`, `TSRHSFunction()`, `DMDATSRHSJacobianLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()` 78714d0ab18SJacob Faibussowitsch S*/ 78814d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *f, void *ctx); 78914d0ab18SJacob Faibussowitsch 79014d0ab18SJacob Faibussowitsch /*S 79114d0ab18SJacob Faibussowitsch DMDATSRHSJacobianLocal - A local residual evaluation function for use with `DMDA` 79214d0ab18SJacob Faibussowitsch 79314d0ab18SJacob Faibussowitsch Calling Sequence: 79414d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 79514d0ab18SJacob Faibussowitsch . t - time at which to evaluate residual 79614d0ab18SJacob Faibussowitsch . x - array of local state information 79714d0ab18SJacob Faibussowitsch . J - Jacobian matrix 79814d0ab18SJacob Faibussowitsch . B - matrix from which to construct the preconditioner; often same as `J` 79914d0ab18SJacob Faibussowitsch - ctx - optional context 80014d0ab18SJacob Faibussowitsch 80114d0ab18SJacob Faibussowitsch Level: beginner 80214d0ab18SJacob Faibussowitsch 80314d0ab18SJacob Faibussowitsch .seealso: `DMDA`, `DMDATSSetRHSJacobianLocal()`, `TSRHSJacobian()`, `DMDATSRHSFunctionLocal()`, `DMDATSIJacobianLocal()`, `DMDATSIFunctionLocal()` 80414d0ab18SJacob Faibussowitsch S*/ 80514d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, Mat J, Mat B, void *ctx); 80614d0ab18SJacob Faibussowitsch 80714d0ab18SJacob Faibussowitsch /*S 80814d0ab18SJacob Faibussowitsch DMDATSIFunctionLocal - A local residual evaluation function for use with `DMDA` 80914d0ab18SJacob Faibussowitsch 81014d0ab18SJacob Faibussowitsch Calling Sequence: 81114d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 81214d0ab18SJacob Faibussowitsch . t - time at which to evaluate residual 81314d0ab18SJacob Faibussowitsch . x - array of local state information 81414d0ab18SJacob Faibussowitsch . xdot - array of local time derivative information 81514d0ab18SJacob Faibussowitsch . imode - output array of local function evaluation information 81614d0ab18SJacob Faibussowitsch - ctx - optional context 81714d0ab18SJacob Faibussowitsch 81814d0ab18SJacob Faibussowitsch Level: beginner 81914d0ab18SJacob Faibussowitsch 82014d0ab18SJacob Faibussowitsch .seealso: `DMDA`, `DMDATSSetIFunctionLocal()`, `DMDATSIJacobianLocal()`, `TSIFunction()` 82114d0ab18SJacob Faibussowitsch S*/ 82214d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, void *imode, void *ctx); 82314d0ab18SJacob Faibussowitsch 82414d0ab18SJacob Faibussowitsch /*S 82514d0ab18SJacob Faibussowitsch DMDATSIJacobianLocal - A local residual evaluation function for use with `DMDA` 82614d0ab18SJacob Faibussowitsch 82714d0ab18SJacob Faibussowitsch Calling Sequence: 82814d0ab18SJacob Faibussowitsch + info - defines the subdomain to evaluate the residual on 82914d0ab18SJacob Faibussowitsch . t - time at which to evaluate the jacobian 83014d0ab18SJacob Faibussowitsch . x - array of local state information 83114d0ab18SJacob Faibussowitsch . xdot - time derivative at this state 83214d0ab18SJacob Faibussowitsch . shift - see `TSSetIJacobian()` for the meaning of this parameter 83314d0ab18SJacob Faibussowitsch . J - Jacobian matrix 83414d0ab18SJacob Faibussowitsch . B - matrix from which to construct the preconditioner; often same as `J` 83514d0ab18SJacob Faibussowitsch - ctx - optional context 83614d0ab18SJacob Faibussowitsch 83714d0ab18SJacob Faibussowitsch Level: beginner 83814d0ab18SJacob Faibussowitsch 83914d0ab18SJacob Faibussowitsch .seealso: `DMDA` `DMDATSSetIJacobianLocal()`, `TSIJacobian()`, `DMDATSIFunctionLocal()`, `DMDATSRHSFunctionLocal()`, `DMDATSRHSJacob ianlocal()` 84014d0ab18SJacob Faibussowitsch S*/ 84114d0ab18SJacob Faibussowitsch PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo *info, PetscReal t, void *x, void *xdot, PetscReal shift, Mat J, Mat B, void *ctx); 84214d0ab18SJacob Faibussowitsch 84314d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM, InsertMode, DMDATSRHSFunctionLocal, void *); 84414d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM, DMDATSRHSJacobianLocal, void *); 84514d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM, InsertMode, DMDATSIFunctionLocal, void *); 84614d0ab18SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM, DMDATSIJacobianLocal, void *); 8476c6b9e74SPeter Brune 848be1b0d75SMatthew G. Knepley typedef struct _n_TSMonitorLGCtx *TSMonitorLGCtx; 849d1212d36SBarry Smith typedef struct { 850d1212d36SBarry Smith Vec ray; 851d1212d36SBarry Smith VecScatter scatter; 852d1212d36SBarry Smith PetscViewer viewer; 85351b4a12fSMatthew G. Knepley TSMonitorLGCtx lgctx; 854d1212d36SBarry Smith } TSMonitorDMDARayCtx; 855d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void **); 856d1212d36SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS, PetscInt, PetscReal, Vec, void *); 85751b4a12fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS, PetscInt, PetscReal, Vec, void *); 858d1212d36SBarry Smith 859bdad233fSMatthew Knepley /* Dynamic creation and loading functions */ 860140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList TSList; 86119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSGetType(TS, TSType *); 86219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSSetType(TS, TSType); 863bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 86430de9b25SBarry Smith 865014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetSNES(TS, SNES *); 866deb2cd25SJed Brown PETSC_EXTERN PetscErrorCode TSSetSNES(TS, SNES); 867014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetKSP(TS, KSP *); 868818ad0c1SBarry Smith 869014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSView(TS, PetscViewer); 87055849f57SBarry Smith PETSC_EXTERN PetscErrorCode TSLoad(TS, PetscViewer); 871fe2efc57SMark PETSC_EXTERN PetscErrorCode TSViewFromOptions(TS, PetscObject, const char[]); 872fe2efc57SMark PETSC_EXTERN PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory, PetscObject, const char[]); 87355849f57SBarry Smith 87455849f57SBarry Smith #define TS_FILE_CLASSID 1211225 87521c89e3eSBarry Smith 876014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS, void *); 877014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS, void *); 87821c89e3eSBarry Smith 879a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtx *); 880a9f9c1f6SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx *); 8814f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS, PetscInt, PetscReal, Vec, void *); 8824f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS, PetscInt, PetscReal, Vec, void *); 88331152f8aSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS, const char *const *); 884b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS, const char *const **); 885e673d494SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx, const char *const *); 886a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS, const char *const *); 887a66092f1SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx, const char *const *); 8887684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 8897684fa3eSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx, PetscErrorCode (*)(void *, Vec, Vec *), PetscErrorCode (*)(void *), void *); 8904f09c107SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS, PetscInt, PetscReal, Vec, void *); 891201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS, PetscInt, PetscReal, Vec, void *); 892201da799SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS, PetscInt, PetscReal, Vec, void *); 893edbaebb3SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorError(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 894d0c080abSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSDMSwarmMonitorMoments(TS, PetscInt, PetscReal, Vec, PetscViewerAndFormat *); 895ef20d060SBarry Smith 896e669de00SBarry Smith struct _n_TSMonitorLGCtxNetwork { 897e669de00SBarry Smith PetscInt nlg; 898e669de00SBarry Smith PetscDrawLG *lg; 899e669de00SBarry Smith PetscBool semilogy; 900e669de00SBarry Smith PetscInt howoften; /* when > 0 uses step % howoften, when negative only final solution plotted */ 901e669de00SBarry Smith }; 902e669de00SBarry Smith typedef struct _n_TSMonitorLGCtxNetwork *TSMonitorLGCtxNetwork; 903e669de00SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *); 904e669de00SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkCreate(TS, const char[], const char[], int, int, int, int, PetscInt, TSMonitorLGCtxNetwork *); 905e669de00SBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorLGCtxNetworkSolution(TS, PetscInt, PetscReal, Vec, void *); 906e669de00SBarry Smith 907b3d3934dSBarry Smith typedef struct _n_TSMonitorEnvelopeCtx *TSMonitorEnvelopeCtx; 908b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS, TSMonitorEnvelopeCtx *); 909b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS, PetscInt, PetscReal, Vec, void *); 910b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS, Vec *, Vec *); 911b3d3934dSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx *); 912b3d3934dSBarry Smith 9138189c53fSBarry Smith typedef struct _n_TSMonitorSPEigCtx *TSMonitorSPEigCtx; 9148189c53fSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, TSMonitorSPEigCtx *); 9158189c53fSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx *); 9168189c53fSBarry Smith PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS, PetscInt, PetscReal, Vec, void *); 9173914022bSBarry Smith 9181b575b74SJoseph Pusztay typedef struct _n_TSMonitorSPCtx *TSMonitorSPCtx; 91960e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscBool, PetscBool, TSMonitorSPCtx *); 9201b575b74SJoseph Pusztay PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx *); 9210ec8ee2bSJoseph Pusztay PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 9221b575b74SJoseph Pusztay 92360e16b1bSMatthew G. Knepley typedef struct _n_TSMonitorHGCtx *TSMonitorHGCtx; 92460e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGCtxCreate(MPI_Comm, const char[], const char[], int, int, int, int, PetscInt, PetscInt, PetscInt, PetscBool, TSMonitorHGCtx *); 92560e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 92660e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGCtxDestroy(TSMonitorHGCtx *); 92760e16b1bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSMonitorHGSwarmSolution(TS, PetscInt, PetscReal, Vec, void *); 92860e16b1bSMatthew G. Knepley 929ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS, PetscInt, PetscInt[], PetscBool[], PetscErrorCode (*)(TS, PetscReal, Vec, PetscReal[], void *), PetscErrorCode (*)(TS, PetscInt, PetscInt[], PetscReal, Vec, PetscBool, void *), void *); 930ca4445c7SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventStep(TS, PetscReal); 931*fe4ad979SIlya Fursov PETSC_EXTERN PetscErrorCode TSSetPostEventSecondStep(TS, PetscReal); 932*fe4ad979SIlya Fursov PETSC_DEPRECATED_FUNCTION(3, 21, 0, "TSSetPostEventSecondStep()", ) static inline PetscErrorCode TSSetPostEventIntervalStep(TS ts, PetscReal dt) 933*fe4ad979SIlya Fursov { 934*fe4ad979SIlya Fursov return TSSetPostEventSecondStep(ts, dt); 935*fe4ad979SIlya Fursov } 9366427ac75SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS, PetscReal, PetscReal[]); 9371ea83e56SMiguel PETSC_EXTERN PetscErrorCode TSGetNumEvents(TS, PetscInt *); 9381ea83e56SMiguel 93976bdecfbSBarry Smith /*J 940195e9b02SBarry Smith TSSSPType - string with the name of a `TSSSP` scheme. 941815f1ad5SJed Brown 942815f1ad5SJed Brown Level: beginner 943815f1ad5SJed Brown 9441cc06b55SBarry Smith .seealso: [](ch_ts), `TSSSPSetType()`, `TS`, `TSSSP` 94576bdecfbSBarry Smith J*/ 94619fd82e9SBarry Smith typedef const char *TSSSPType; 947815f1ad5SJed Brown #define TSSSPRKS2 "rks2" 948815f1ad5SJed Brown #define TSSSPRKS3 "rks3" 949815f1ad5SJed Brown #define TSSSPRK104 "rk104" 950815f1ad5SJed Brown 95119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSSSPSetType(TS, TSSSPType); 95219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSSSPGetType(TS, TSSSPType *); 953014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS, PetscInt); 954014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS, PetscInt *); 955787849ffSJed Brown PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 956560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 957013797aaSBarry Smith PETSC_EXTERN PetscFunctionList TSSSPList; 958815f1ad5SJed Brown 959ed657a08SJed Brown /*S 96084df9cb4SJed Brown TSAdapt - Abstract object that manages time-step adaptivity 96184df9cb4SJed Brown 96284df9cb4SJed Brown Level: beginner 96384df9cb4SJed Brown 9641cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdaptCreate()`, `TSAdaptType` 96584df9cb4SJed Brown S*/ 96684df9cb4SJed Brown typedef struct _p_TSAdapt *TSAdapt; 96784df9cb4SJed Brown 968f0fc11ceSJed Brown /*J 96987497f52SBarry Smith TSAdaptType - String with the name of `TSAdapt` scheme. 97084df9cb4SJed Brown 97184df9cb4SJed Brown Level: beginner 97284df9cb4SJed Brown 9731cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptSetType()`, `TS`, `TSAdapt` 974f0fc11ceSJed Brown J*/ 975f8963224SJed Brown typedef const char *TSAdaptType; 976b0f836d7SJed Brown #define TSADAPTNONE "none" 9771917a363SLisandro Dalcin #define TSADAPTBASIC "basic" 978cb7ab186SLisandro Dalcin #define TSADAPTDSP "dsp" 9798d59e960SJed Brown #define TSADAPTCFL "cfl" 9801917a363SLisandro Dalcin #define TSADAPTGLEE "glee" 981d949e4a4SStefano Zampini #define TSADAPTHISTORY "history" 98284df9cb4SJed Brown 983552698daSJed Brown PETSC_EXTERN PetscErrorCode TSGetAdapt(TS, TSAdapt *); 984bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[], PetscErrorCode (*)(TSAdapt)); 985607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 986014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 987014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm, TSAdapt *); 98819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt, TSAdaptType); 989d0288e4fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt, TSAdaptType *); 990014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt, const char[]); 991014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 992014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt, const char[], PetscInt, PetscInt, PetscReal, PetscReal, PetscBool); 993014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt, PetscInt *, const PetscInt **, const PetscInt **, const PetscReal **, const PetscReal **); 994014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt, TS, PetscReal, PetscInt *, PetscReal *, PetscBool *); 995b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt, TS, PetscReal, Vec, PetscBool *); 996014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt, PetscViewer); 997ad6bc421SBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt, PetscViewer); 998dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt, PetscOptionItems *); 999099af0a3SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt); 1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt *); 1001014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt, PetscBool); 1002bf997491SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt, PetscBool); 10031917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt, PetscReal, PetscReal); 10041917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt, PetscReal *, PetscReal *); 100576cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt, PetscReal); 100676cddca1SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt, PetscReal *); 10071917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt, PetscReal, PetscReal); 10081917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt, PetscReal *, PetscReal *); 100962c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt, PetscReal); 101062c23b28SMark PETSC_EXTERN PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt, PetscReal *); 1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt, PetscReal, PetscReal); 10121917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt, PetscReal *, PetscReal *); 1013b295832fSPierre Barbier de Reuille PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt, PetscErrorCode (*)(TSAdapt, TS, PetscReal, Vec, PetscBool *)); 1014d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt, PetscInt n, PetscReal hist[], PetscBool); 101575017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt, TSTrajectory, PetscBool); 101675017583SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt, PetscInt, PetscReal *, PetscReal *); 1017de50f1caSBarry Smith PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt, PetscInt); 1018cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt, const char *); 1019cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt, PetscReal, PetscReal, PetscReal); 102084df9cb4SJed Brown 102184df9cb4SJed Brown /*S 102287497f52SBarry Smith TSGLLEAdapt - Abstract object that manages time-step adaptivity for `TSGLLE` 1023ed657a08SJed Brown 1024ed657a08SJed Brown Level: beginner 1025ed657a08SJed Brown 1026195e9b02SBarry Smith Developer Note: 102787497f52SBarry Smith This functionality should be replaced by the `TSAdapt`. 102884df9cb4SJed Brown 10291cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdaptCreate()`, `TSGLLEAdaptType` 1030ed657a08SJed Brown S*/ 103126d28e4eSEmil Constantinescu typedef struct _p_TSGLLEAdapt *TSGLLEAdapt; 1032ed657a08SJed Brown 103376bdecfbSBarry Smith /*J 103487497f52SBarry Smith TSGLLEAdaptType - String with the name of `TSGLLEAdapt` scheme 1035ed657a08SJed Brown 1036ed657a08SJed Brown Level: beginner 1037ed657a08SJed Brown 1038195e9b02SBarry Smith Developer Note: 1039195e9b02SBarry Smith This functionality should be replaced by the `TSAdaptType`. 1040195e9b02SBarry Smith 10411cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLEAdaptSetType()`, `TS` 104276bdecfbSBarry Smith J*/ 104326d28e4eSEmil Constantinescu typedef const char *TSGLLEAdaptType; 104426d28e4eSEmil Constantinescu #define TSGLLEADAPT_NONE "none" 104526d28e4eSEmil Constantinescu #define TSGLLEADAPT_SIZE "size" 104626d28e4eSEmil Constantinescu #define TSGLLEADAPT_BOTH "both" 1047ed657a08SJed Brown 104826d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[], PetscErrorCode (*)(TSGLLEAdapt)); 104926d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void); 105026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void); 105126d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm, TSGLLEAdapt *); 105226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt, TSGLLEAdaptType); 105326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt, const char[]); 105426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt, PetscInt, const PetscInt[], const PetscReal[], const PetscReal[], PetscInt, PetscReal, PetscReal, PetscInt *, PetscReal *, PetscBool *); 105526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt, PetscViewer); 1056dbbe0bcdSBarry Smith PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(TSGLLEAdapt, PetscOptionItems *); 105726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt *); 1058ed657a08SJed Brown 105976bdecfbSBarry Smith /*J 106087497f52SBarry Smith TSGLLEAcceptType - String with the name of `TSGLLEAccept` scheme 1061ed657a08SJed Brown 1062ed657a08SJed Brown Level: beginner 1063ed657a08SJed Brown 10641cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLLESetAcceptType()`, `TS`, `TSGLLEAccept` 106576bdecfbSBarry Smith J*/ 106626d28e4eSEmil Constantinescu typedef const char *TSGLLEAcceptType; 106726d28e4eSEmil Constantinescu #define TSGLLEACCEPT_ALWAYS "always" 1068ed657a08SJed Brown 106926d28e4eSEmil Constantinescu PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS, PetscReal, PetscReal, const PetscReal[], PetscBool *); 107026d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[], TSGLLEAcceptFunction); 1071ed657a08SJed Brown 107276bdecfbSBarry Smith /*J 1073195e9b02SBarry Smith TSGLLEType - string with the name of a General Linear `TSGLLE` type 107418b56cb9SJed Brown 107518b56cb9SJed Brown Level: beginner 107618b56cb9SJed Brown 10771cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLESetType()`, `TSGLLERegister()`, `TSGLLEAccept` 107876bdecfbSBarry Smith J*/ 107926d28e4eSEmil Constantinescu typedef const char *TSGLLEType; 108026d28e4eSEmil Constantinescu #define TSGLLE_IRKS "irks" 108118b56cb9SJed Brown 108226d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[], PetscErrorCode (*)(TS)); 108326d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void); 108426d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void); 108526d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetType(TS, TSGLLEType); 108626d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS, TSGLLEAdapt *); 108726d28e4eSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS, TSGLLEAcceptType); 108818b56cb9SJed Brown 1089020d8f30SJed Brown /*J 1090195e9b02SBarry Smith TSEIMEXType - String with the name of an Extrapolated IMEX `TSEIMEX` type 10917d5697caSHong Zhang 10927d5697caSHong Zhang Level: beginner 10937d5697caSHong Zhang 10941cc06b55SBarry Smith .seealso: [](ch_ts), `TSEIMEXSetType()`, `TS`, `TSEIMEX`, `TSEIMEXRegister()` 10957d5697caSHong Zhang J*/ 10967d5697caSHong Zhang #define TSEIMEXType char * 10977d5697caSHong Zhang 10987d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt); 10997d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt, PetscInt); 11007d5697caSHong Zhang PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS, PetscBool); 11017d5697caSHong Zhang 11027d5697caSHong Zhang /*J 1103195e9b02SBarry Smith TSRKType - String with the name of a Runge-Kutta `TSRK` type 1104f68a32c8SEmil Constantinescu 1105f68a32c8SEmil Constantinescu Level: beginner 1106f68a32c8SEmil Constantinescu 11071cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSRKSetType()`, `TS`, `TSRK`, `TSRKRegister()` 1108f68a32c8SEmil Constantinescu J*/ 1109f68a32c8SEmil Constantinescu typedef const char *TSRKType; 1110f68a32c8SEmil Constantinescu #define TSRK1FE "1fe" 1111fdefd5e5SDebojyoti Ghosh #define TSRK2A "2a" 1112e7685601SHong Zhang #define TSRK2B "2b" 1113f68a32c8SEmil Constantinescu #define TSRK3 "3" 1114fdefd5e5SDebojyoti Ghosh #define TSRK3BS "3bs" 1115f68a32c8SEmil Constantinescu #define TSRK4 "4" 1116fdefd5e5SDebojyoti Ghosh #define TSRK5F "5f" 1117fdefd5e5SDebojyoti Ghosh #define TSRK5DP "5dp" 111805e23783SLisandro Dalcin #define TSRK5BS "5bs" 111937acaa02SHendrik Ranocha #define TSRK6VR "6vr" 112037acaa02SHendrik Ranocha #define TSRK7VR "7vr" 112137acaa02SHendrik Ranocha #define TSRK8VR "8vr" 112205e23783SLisandro Dalcin 11232ea87ba9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKGetOrder(TS, PetscInt *); 11240f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetType(TS, TSRKType *); 11250f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKSetType(TS, TSRKType); 11260f7a28cdSStefano Zampini PETSC_EXTERN PetscErrorCode TSRKGetTableau(TS, PetscInt *, const PetscReal **, const PetscReal **, const PetscReal **, const PetscReal **, PetscInt *, const PetscReal **, PetscBool *); 11270fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS, PetscBool); 11280fe4d17eSHong Zhang PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS, PetscBool *); 1129f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1130f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 1131560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 1132f68a32c8SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 1133f68a32c8SEmil Constantinescu 1134f68a32c8SEmil Constantinescu /*J 1135195e9b02SBarry Smith TSMPRKType - String with the name of a Partitioned Runge-Kutta `TSMPRK` type 11369f537349Sluzhanghpp 11379f537349Sluzhanghpp Level: beginner 11389f537349Sluzhanghpp 11391cc06b55SBarry Smith .seealso: [](ch_ts), `TSMPRKSetType()`, `TS`, `TSMPRK`, `TSMPRKRegister()` 11409f537349Sluzhanghpp J*/ 11414b84eec9SHong Zhang typedef const char *TSMPRKType; 114219c2959aSHong Zhang #define TSMPRK2A22 "2a22" 114319c2959aSHong Zhang #define TSMPRK2A23 "2a23" 114419c2959aSHong Zhang #define TSMPRK2A32 "2a32" 114519c2959aSHong Zhang #define TSMPRK2A33 "2a33" 11464b84eec9SHong Zhang #define TSMPRKP2 "p2" 11474b84eec9SHong Zhang #define TSMPRKP3 "p3" 11489f537349Sluzhanghpp 11494b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *); 11504b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType); 115179bc8a77SHong 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[]); 11524b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void); 11534b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void); 11544b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void); 11559f537349Sluzhanghpp 11569f537349Sluzhanghpp /*J 1157195e9b02SBarry Smith TSIRKType - String with the name of an implicit Runge-Kutta `TSIRK` type 1158d2567f34SHong Zhang 1159d2567f34SHong Zhang Level: beginner 1160d2567f34SHong Zhang 11611cc06b55SBarry Smith .seealso: [](ch_ts), `TSIRKSetType()`, `TS`, `TSIRK`, `TSIRKRegister()` 1162d2567f34SHong Zhang J*/ 1163d2567f34SHong Zhang typedef const char *TSIRKType; 1164d2567f34SHong Zhang #define TSIRKGAUSS "gauss" 1165d2567f34SHong Zhang 1166d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetType(TS, TSIRKType *); 1167d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetType(TS, TSIRKType); 1168d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKGetNumStages(TS, PetscInt *); 1169d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKSetNumStages(TS, PetscInt); 1170d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegister(const char[], PetscErrorCode (*function)(TS)); 1171d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKTableauCreate(TS, PetscInt, const PetscReal *, const PetscReal *, const PetscReal *, const PetscReal *, const PetscScalar *, const PetscScalar *, const PetscScalar *); 1172d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKInitializePackage(void); 1173d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKFinalizePackage(void); 1174d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSIRKRegisterDestroy(void); 1175d2567f34SHong Zhang 1176d2567f34SHong Zhang /*J 1177195e9b02SBarry Smith TSGLEEType - String with the name of a General Linear with Error Estimation `TSGLEE` type 1178b6a60446SDebojyoti Ghosh 1179b6a60446SDebojyoti Ghosh Level: beginner 1180b6a60446SDebojyoti Ghosh 11811cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetType()`, `TS`, `TSGLEE`, `TSGLEERegister()` 1182b6a60446SDebojyoti Ghosh J*/ 1183b6a60446SDebojyoti Ghosh typedef const char *TSGLEEType; 1184ab8c5912SEmil Constantinescu #define TSGLEEi1 "BE1" 1185b6a60446SDebojyoti Ghosh #define TSGLEE23 "23" 1186b6a60446SDebojyoti Ghosh #define TSGLEE24 "24" 1187b6a60446SDebojyoti Ghosh #define TSGLEE25I "25i" 1188b6a60446SDebojyoti Ghosh #define TSGLEE35 "35" 1189b6a60446SDebojyoti Ghosh #define TSGLEEEXRK2A "exrk2a" 1190b6a60446SDebojyoti Ghosh #define TSGLEERK32G1 "rk32g1" 1191b6a60446SDebojyoti Ghosh #define TSGLEERK285EX "rk285ex" 119216a05f60SBarry Smith 1193b6a60446SDebojyoti Ghosh /*J 1194195e9b02SBarry Smith TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation `TSGLEE` type 1195b6a60446SDebojyoti Ghosh 1196b6a60446SDebojyoti Ghosh Level: beginner 1197b6a60446SDebojyoti Ghosh 11981cc06b55SBarry Smith .seealso: [](ch_ts), `TSGLEESetMode()`, `TS`, `TSGLEE`, `TSGLEERegister()` 1199b6a60446SDebojyoti Ghosh J*/ 1200b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts, TSGLEEType *); 1201b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts, TSGLEEType); 120257df6a1bSDebojyoti 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[]); 12034bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSGLEERegisterAll(void); 1204b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void); 1205b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void); 1206b6a60446SDebojyoti Ghosh PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void); 1207b6a60446SDebojyoti Ghosh 1208b6a60446SDebojyoti Ghosh /*J 1209195e9b02SBarry Smith TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX `TSARKIMEX` type 1210020d8f30SJed Brown 1211020d8f30SJed Brown Level: beginner 1212020d8f30SJed Brown 12131cc06b55SBarry Smith .seealso: [](ch_ts), `TSARKIMEXSetType()`, `TS`, `TSARKIMEX`, `TSARKIMEXRegister()` 1214020d8f30SJed Brown J*/ 121519fd82e9SBarry Smith typedef const char *TSARKIMEXType; 1216e817cc15SEmil Constantinescu #define TSARKIMEX1BEE "1bee" 12171f80e275SEmil Constantinescu #define TSARKIMEXA2 "a2" 12181f80e275SEmil Constantinescu #define TSARKIMEXL2 "l2" 12191f80e275SEmil Constantinescu #define TSARKIMEXARS122 "ars122" 12201f80e275SEmil Constantinescu #define TSARKIMEX2C "2c" 12218a381b04SJed Brown #define TSARKIMEX2D "2d" 1222a3a57f36SJed Brown #define TSARKIMEX2E "2e" 12236cf0794eSJed Brown #define TSARKIMEXPRSSP2 "prssp2" 12248a381b04SJed Brown #define TSARKIMEX3 "3" 12256cf0794eSJed Brown #define TSARKIMEXBPR3 "bpr3" 12266cf0794eSJed Brown #define TSARKIMEXARS443 "ars443" 12278a381b04SJed Brown #define TSARKIMEX4 "4" 12288a381b04SJed Brown #define TSARKIMEX5 "5" 122919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts, TSARKIMEXType *); 123019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts, TSARKIMEXType); 1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS, PetscBool); 12323a28c0e5SStefano Zampini PETSC_EXTERN PetscErrorCode TSARKIMEXGetFullyImplicit(TS, PetscBool *); 123319fd82e9SBarry 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[]); 1234607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 1235560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 12378a381b04SJed Brown 1238020d8f30SJed Brown /*J 1239d27334e2SStefano Zampini TSDIRKType - String with the name of a Diagonally Implicit Runge-Kutta `TSDIRK` type 1240d27334e2SStefano Zampini 1241d27334e2SStefano Zampini Level: beginner 1242d27334e2SStefano Zampini 1243d27334e2SStefano Zampini .seealso: [](ch_ts), `TSDIRKSetType()`, `TS`, `TSDIRK`, `TSDIRKRegister()` 1244d27334e2SStefano Zampini J*/ 1245d27334e2SStefano Zampini typedef const char *TSDIRKType; 1246d27334e2SStefano Zampini #define TSDIRKS212 "s212" 12473405e92cSStefano Zampini #define TSDIRKES122SAL "es122sal" 1248d27334e2SStefano Zampini #define TSDIRKES213SAL "es213sal" 1249d27334e2SStefano Zampini #define TSDIRKES324SAL "es324sal" 1250d27334e2SStefano Zampini #define TSDIRKES325SAL "es325sal" 1251d27334e2SStefano Zampini #define TSDIRK657A "657a" 1252d27334e2SStefano Zampini #define TSDIRKES648SA "es648sa" 1253d27334e2SStefano Zampini #define TSDIRK658A "658a" 1254d27334e2SStefano Zampini #define TSDIRKS659A "s659a" 1255d27334e2SStefano Zampini #define TSDIRK7510SAL "7510sal" 1256d27334e2SStefano Zampini #define TSDIRKES7510SA "es7510sa" 1257d27334e2SStefano Zampini #define TSDIRK759A "759a" 1258d27334e2SStefano Zampini #define TSDIRKS7511SAL "s7511sal" 1259d27334e2SStefano Zampini #define TSDIRK8614A "8614a" 1260d27334e2SStefano Zampini #define TSDIRK8616SAL "8616sal" 1261d27334e2SStefano Zampini #define TSDIRKES8516SAL "es8516sal" 1262d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKGetType(TS ts, TSDIRKType *); 1263d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKSetType(TS ts, TSDIRKType); 1264d27334e2SStefano Zampini PETSC_EXTERN PetscErrorCode TSDIRKRegister(TSDIRKType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 1265d27334e2SStefano Zampini 1266d27334e2SStefano Zampini /*J 1267195e9b02SBarry Smith TSRosWType - String with the name of a Rosenbrock-W `TSROSW` type 1268020d8f30SJed Brown 1269020d8f30SJed Brown Level: beginner 1270020d8f30SJed Brown 12711cc06b55SBarry Smith .seealso: [](ch_ts), `TSRosWSetType()`, `TS`, `TSROSW`, `TSRosWRegister()` 1272020d8f30SJed Brown J*/ 127319fd82e9SBarry Smith typedef const char *TSRosWType; 127461692a83SJed Brown #define TSROSW2M "2m" 127561692a83SJed Brown #define TSROSW2P "2p" 1276fe7e6d57SJed Brown #define TSROSWRA3PW "ra3pw" 1277fe7e6d57SJed Brown #define TSROSWRA34PW2 "ra34pw2" 1278ef3c5b88SJed Brown #define TSROSWRODAS3 "rodas3" 1279ef3c5b88SJed Brown #define TSROSWSANDU3 "sandu3" 1280961f28d0SJed Brown #define TSROSWASSP3P3S1C "assp3p3s1c" 1281961f28d0SJed Brown #define TSROSWLASSP3P4S2C "lassp3p4s2c" 128243b21953SEmil Constantinescu #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 1283753f8adbSEmil Constantinescu #define TSROSWARK3 "ark3" 12843606a31eSEmil Constantinescu #define TSROSWTHETA1 "theta1" 12853606a31eSEmil Constantinescu #define TSROSWTHETA2 "theta2" 128642faf41dSJed Brown #define TSROSWGRK4T "grk4t" 128742faf41dSJed Brown #define TSROSWSHAMP4 "shamp4" 128842faf41dSJed Brown #define TSROSWVELDD4 "veldd4" 128942faf41dSJed Brown #define TSROSW4L "4l" 1290b1c69cc3SEmil Constantinescu 12916bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWGetType(TS, TSRosWType *); 12926bd7aeb5SHong Zhang PETSC_EXTERN PetscErrorCode TSRosWSetType(TS, TSRosWType); 1293014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS, PetscBool); 129419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType, PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], PetscInt, const PetscReal[]); 129519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal); 1296607a6623SBarry Smith PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 1297560360afSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 1298014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 1299e27a552bSJed Brown 1300211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS, PetscInt); 1301211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS, PetscInt *); 1302211a84d6SLisandro Dalcin 13036bd7aeb5SHong Zhang /*J 1304195e9b02SBarry Smith TSBasicSymplecticType - String with the name of a basic symplectic integration `TSBASICSYMPLECTIC` type 13056bd7aeb5SHong Zhang 13066bd7aeb5SHong Zhang Level: beginner 13076bd7aeb5SHong Zhang 13081cc06b55SBarry Smith .seealso: [](ch_ts), `TSBasicSymplecticSetType()`, `TS`, `TSBASICSYMPLECTIC`, `TSBasicSymplecticRegister()` 13096bd7aeb5SHong Zhang J*/ 1310e49d4f37SHong Zhang typedef const char *TSBasicSymplecticType; 1311e49d4f37SHong Zhang #define TSBASICSYMPLECTICSIEULER "1" 1312e49d4f37SHong Zhang #define TSBASICSYMPLECTICVELVERLET "2" 1313e49d4f37SHong Zhang #define TSBASICSYMPLECTIC3 "3" 1314e49d4f37SHong Zhang #define TSBASICSYMPLECTIC4 "4" 1315e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS, TSBasicSymplecticType); 1316e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS, TSBasicSymplecticType *); 1317e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType, PetscInt, PetscInt, PetscReal[], PetscReal[]); 13184bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterAll(void); 1319e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void); 1320e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void); 1321e49d4f37SHong Zhang PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void); 13226bd7aeb5SHong Zhang 1323db4653c2SDaniel Finn /*J 1324195e9b02SBarry Smith TSDISCGRAD - The Discrete Gradient integrator is a timestepper for Hamiltonian systems designed to conserve the first integral (energy), 132587497f52SBarry Smith but also has the property for some systems of monotonicity in a functional. 1326db4653c2SDaniel Finn 1327db4653c2SDaniel Finn Level: beginner 1328db4653c2SDaniel Finn 13291cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, TSDiscGradSetFormulation()`, `TSDiscGradGetFormulation()` 1330db4653c2SDaniel Finn J*/ 133140be0ff1SMatthew 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 *); 13324bf303faSJacob 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 *); 1333db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradIsGonzalez(TS, PetscBool *); 1334db4653c2SDaniel Finn PETSC_EXTERN PetscErrorCode TSDiscGradUseGonzalez(TS, PetscBool); 133540be0ff1SMatthew G. Knepley 133683e2fdc7SBarry Smith /* 13370f3b3ca1SHong Zhang PETSc interface to Sundials 133883e2fdc7SBarry Smith */ 1339e808b789SPatrick Sanan #ifdef PETSC_HAVE_SUNDIALS2 13409371c9d4SSatish Balay typedef enum { 13419371c9d4SSatish Balay SUNDIALS_ADAMS = 1, 13429371c9d4SSatish Balay SUNDIALS_BDF = 2 13439371c9d4SSatish Balay } TSSundialsLmmType; 13446a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 13459371c9d4SSatish Balay typedef enum { 13469371c9d4SSatish Balay SUNDIALS_MODIFIED_GS = 1, 13479371c9d4SSatish Balay SUNDIALS_CLASSICAL_GS = 2 13489371c9d4SSatish Balay } TSSundialsGramSchmidtType; 13496a6fc655SJed Brown PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 13504bf303faSJacob Faibussowitsch 1351014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS, TSSundialsLmmType); 1352014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS, PC *); 1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS, PetscReal, PetscReal); 1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS, PetscReal); 1355014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS, PetscReal); 1356014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS, PetscInt *, PetscInt *); 1357014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS, TSSundialsGramSchmidtType); 1358014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS, PetscInt); 1359014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS, PetscReal); 1360014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS, PetscBool); 1361014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS, PetscInt); 1362c4e80e11SFlorian PETSC_EXTERN PetscErrorCode TSSundialsSetMaxord(TS, PetscInt); 1363918687b8SPatrick Sanan PETSC_EXTERN PetscErrorCode TSSundialsSetUseDense(TS, PetscBool); 13646dc17ff5SMatthew Knepley #endif 136583e2fdc7SBarry Smith 1366014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS, PetscReal); 1367014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS, PetscReal *); 1368014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS, PetscBool *); 1369014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS, PetscBool); 13700de4c49aSJed Brown 1371014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS, PetscReal); 1372014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS, PetscReal, PetscReal, PetscReal); 1373014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS, PetscReal *, PetscReal *, PetscReal *); 137488df8a41SLisandro Dalcin 1375818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS, PetscReal); 1376818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS, PetscReal, PetscReal, PetscReal, PetscReal); 1377818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 1378818efac9SLisandro Dalcin 1379220f924aSDavid Kamensky /*S 1380220f924aSDavid Kamensky TSAlpha2Predictor - A callback to set the predictor (i.e., the initial guess for the nonlinear solver) in 1381220f924aSDavid Kamensky a second-order generalized-alpha time integrator. 1382220f924aSDavid Kamensky 1383220f924aSDavid Kamensky Calling Sequence: 1384220f924aSDavid Kamensky + ts - the `TS` context obtained from `TSCreate()` 1385220f924aSDavid Kamensky . X0 - the previous time step's state vector 1386220f924aSDavid Kamensky . V0 - the previous time step's first derivative of the state vector 1387220f924aSDavid Kamensky . A0 - the previous time step's second derivative of the state vector 1388220f924aSDavid Kamensky . X1 - the vector into which the initial guess for the current time step will be written 1389220f924aSDavid Kamensky - ctx - [optional] user-defined context for the predictor evaluation routine (may be `NULL`) 1390220f924aSDavid Kamensky 1391220f924aSDavid Kamensky Level: intermediate 1392220f924aSDavid Kamensky 1393220f924aSDavid Kamensky .seealso: [](ch_ts), `TS`, `TSAlpha2SetPredictor()` 1394220f924aSDavid Kamensky S*/ 1395220f924aSDavid Kamensky PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSAlpha2Predictor)(TS ts, Vec X0, Vec V0, Vec A0, Vec X1, void *ctx); 1396220f924aSDavid Kamensky PETSC_EXTERN PetscErrorCode TSAlpha2SetPredictor(TS, TSAlpha2Predictor, void *ctx); 1397220f924aSDavid Kamensky 1398014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSSetDM(TS, DM); 1399014dd563SJed Brown PETSC_EXTERN PetscErrorCode TSGetDM(TS, DM *); 14006c699258SBarry Smith 1401014dd563SJed Brown PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES, Vec, Vec, void *); 1402d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES, Vec, Mat, Mat, void *); 14030f5c6efeSJed Brown 1404f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS, PetscBool *); 1405f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS, PetscBool *); 1406aad739acSMatthew G. Knepley 14072e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeInitialCondition(TS, PetscErrorCode (**)(TS, Vec)); 14082e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeInitialCondition(TS, PetscErrorCode (*)(TS, Vec)); 14092e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeInitialCondition(TS, Vec); 14102e61be88SMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSGetComputeExactError(TS, PetscErrorCode (**)(TS, Vec, Vec)); 1411aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSSetComputeExactError(TS, PetscErrorCode (*)(TS, Vec, Vec)); 1412aad739acSMatthew G. Knepley PETSC_EXTERN PetscErrorCode TSComputeExactError(TS, Vec, Vec); 1413f2ed2dc7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscConvEstUseTS(PetscConvEst, PetscBool); 1414d60b7d5cSBarry Smith 1415d60b7d5cSBarry Smith PETSC_EXTERN PetscErrorCode TSSetMatStructure(TS, MatStructure); 1416