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