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