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 TSSetCostHessianProducts(TS,PetscInt,Vec*,Vec*,Vec); 273 PETSC_EXTERN PetscErrorCode TSGetCostHessianProducts(TS,PetscInt*,Vec**,Vec**,Vec*); 274 275 276 /*S 277 TSTrajectory - Abstract PETSc object that stores the trajectory (solution of ODE/DAE at each time step) 278 279 Level: advanced 280 281 Concepts: ODE solvers, trajectory 282 283 .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectorySetType(), TSTrajectoryDestroy(), TSTrajectoryReset() 284 S*/ 285 typedef struct _p_TSTrajectory* TSTrajectory; 286 287 /*J 288 TSTrajectorySetType - String with the name of a PETSc TS trajectory storage method 289 290 Level: intermediate 291 292 .seealso: TSSetSaveTrajectory(), TSTrajectoryCreate(), TSTrajectoryDestroy() 293 J*/ 294 typedef const char* TSTrajectoryType; 295 #define TSTRAJECTORYBASIC "basic" 296 #define TSTRAJECTORYSINGLEFILE "singlefile" 297 #define TSTRAJECTORYMEMORY "memory" 298 #define TSTRAJECTORYVISUALIZATION "visualization" 299 300 PETSC_EXTERN PetscFunctionList TSTrajectoryList; 301 PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID; 302 PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled; 303 304 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS); 305 PETSC_EXTERN PetscErrorCode TSResetTrajectory(TS); 306 307 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*); 308 PETSC_EXTERN PetscErrorCode TSTrajectoryReset(TSTrajectory); 309 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*); 310 PETSC_EXTERN PetscErrorCode TSTrajectoryView(TSTrajectory,PetscViewer); 311 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,TSTrajectoryType); 312 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec); 313 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*); 314 PETSC_EXTERN PetscErrorCode TSTrajectoryGetVecs(TSTrajectory,TS,PetscInt,PetscReal*,Vec,Vec); 315 PETSC_EXTERN PetscErrorCode TSTrajectoryGetUpdatedHistoryVecs(TSTrajectory,TS,PetscReal,Vec*,Vec*); 316 PETSC_EXTERN PetscErrorCode TSTrajectoryGetNumSteps(TSTrajectory,PetscInt*); 317 PETSC_EXTERN PetscErrorCode TSTrajectoryRestoreUpdatedHistoryVecs(TSTrajectory,Vec*,Vec*); 318 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS); 319 PETSC_EXTERN PetscErrorCode TSTrajectoryRegister(const char[],PetscErrorCode (*)(TSTrajectory,TS)); 320 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void); 321 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS); 322 PETSC_EXTERN PetscErrorCode TSTrajectorySetMonitor(TSTrajectory,PetscBool); 323 PETSC_EXTERN PetscErrorCode TSTrajectorySetVariableNames(TSTrajectory,const char * const*); 324 PETSC_EXTERN PetscErrorCode TSTrajectorySetTransform(TSTrajectory,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*); 325 PETSC_EXTERN PetscErrorCode TSTrajectorySetSolutionOnly(TSTrajectory,PetscBool); 326 PETSC_EXTERN PetscErrorCode TSTrajectoryGetSolutionOnly(TSTrajectory,PetscBool*); 327 PETSC_EXTERN PetscErrorCode TSTrajectorySetKeepFiles(TSTrajectory,PetscBool); 328 PETSC_EXTERN PetscErrorCode TSTrajectorySetDirname(TSTrajectory,const char[]); 329 PETSC_EXTERN PetscErrorCode TSTrajectorySetFiletemplate(TSTrajectory,const char[]); 330 PETSC_EXTERN PetscErrorCode TSGetTrajectory(TS,TSTrajectory*); 331 332 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*); 333 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**); 334 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*); 335 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*); 336 PETSC_EXTERN PetscErrorCode TSComputeCostIntegrand(TS,PetscReal,Vec,Vec); 337 338 PETSC_EXTERN PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems*,TS); 339 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*); 340 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**)); 341 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS); 342 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*)); 343 344 PETSC_EXTERN PETSC_DEPRECATED("Use TSSetRHSJacobianP") PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*); 345 PETSC_EXTERN PETSC_DEPRECATED("Use TSComputeRHSJacobianP") PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat); 346 PETSC_EXTERN PETSC_DEPRECATED("Use TSComputeDRDPFunction") PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*); 347 PETSC_EXTERN PETSC_DEPRECATED("Use TSComputeDRDUFunction") PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*); 348 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS); 349 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt); 350 351 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 352 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 353 PETSC_EXTERN PetscErrorCode TSAdjointReset(TS); 354 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS); 355 PETSC_EXTERN PetscErrorCode TSAdjointInitializeForward(TS,Mat); 356 357 PETSC_EXTERN PetscErrorCode TSForwardSetSensitivities(TS,PetscInt,Mat); 358 PETSC_EXTERN PetscErrorCode TSForwardGetSensitivities(TS,PetscInt*,Mat*); 359 PETSC_EXTERN PetscErrorCode TSForwardSetIntegralGradients(TS,PetscInt,Vec *); 360 PETSC_EXTERN PetscErrorCode TSForwardGetIntegralGradients(TS,PetscInt*,Vec **); 361 PETSC_EXTERN PetscErrorCode TSForwardSetRHSJacobianP(TS,Vec*,PetscErrorCode(*)(TS,PetscReal,Vec,Vec*,void*),void*); 362 PETSC_EXTERN PetscErrorCode TSForwardComputeRHSJacobianP(TS,PetscReal,Vec,Vec*); 363 PETSC_EXTERN PetscErrorCode TSForwardSetUp(TS); 364 PETSC_EXTERN PetscErrorCode TSForwardReset(TS); 365 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS); 366 PETSC_EXTERN PetscErrorCode TSForwardStep(TS); 367 PETSC_EXTERN PetscErrorCode TSForwardSetInitialSensitivities(TS,Mat); 368 PETSC_EXTERN PetscErrorCode TSForwardGetStages(TS,PetscInt*,Mat**); 369 370 PETSC_EXTERN PetscErrorCode TSSetMaxSteps(TS,PetscInt); 371 PETSC_EXTERN PetscErrorCode TSGetMaxSteps(TS,PetscInt*); 372 PETSC_EXTERN PetscErrorCode TSSetMaxTime(TS,PetscReal); 373 PETSC_EXTERN PetscErrorCode TSGetMaxTime(TS,PetscReal*); 374 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 375 PETSC_EXTERN PetscErrorCode TSGetExactFinalTime(TS,TSExactFinalTimeOption*); 376 377 PETSC_EXTERN PETSC_DEPRECATED("Use TSSetTime[Step]") PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 378 PETSC_EXTERN PETSC_DEPRECATED("Use TSSetMax{Steps|Time}") PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 379 PETSC_EXTERN PETSC_DEPRECATED("Use TSGetMax{Steps|Time}") PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 380 PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber") PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 381 PETSC_EXTERN PETSC_DEPRECATED("Use TSGetStepNumber") PetscErrorCode TSGetTotalSteps(TS,PetscInt*); 382 383 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*); 384 PETSC_EXTERN PetscErrorCode TSMonitorExtreme(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*); 385 386 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 387 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 388 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 389 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 390 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*); 391 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 392 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionFunction(TS,PetscInt,PetscReal,Vec,void*); 393 394 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *); 395 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*); 396 397 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*); 398 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 399 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 400 401 PETSC_EXTERN PetscErrorCode TSStep(TS); 402 PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*); 403 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 404 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 405 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*); 406 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType); 407 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 408 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 409 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 410 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 411 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 412 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 413 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 414 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 415 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 416 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 417 PETSC_EXTERN PetscErrorCode TSRestartStep(TS); 418 PETSC_EXTERN PetscErrorCode TSRollBack(TS); 419 420 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**); 421 422 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 423 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 424 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*); 425 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 426 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 427 PETSC_EXTERN PetscErrorCode TSGetStepNumber(TS,PetscInt*); 428 PETSC_EXTERN PetscErrorCode TSSetStepNumber(TS,PetscInt); 429 430 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 431 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*); 432 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 433 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 434 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 435 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 436 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool); 437 438 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 439 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 440 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSForcingFunction)(TS,PetscReal,Vec,void*); 441 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,TSForcingFunction,void*); 442 443 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 444 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 445 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 446 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 447 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 448 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 449 450 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Function)(TS,PetscReal,Vec,Vec,Vec,Vec,void*); 451 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSI2Jacobian)(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat,void*); 452 PETSC_EXTERN PetscErrorCode TSSetI2Function(TS,Vec,TSI2Function,void*); 453 PETSC_EXTERN PetscErrorCode TSGetI2Function(TS,Vec*,TSI2Function*,void**); 454 PETSC_EXTERN PetscErrorCode TSSetI2Jacobian(TS,Mat,Mat,TSI2Jacobian,void*); 455 PETSC_EXTERN PetscErrorCode TSGetI2Jacobian(TS,Mat*,Mat*,TSI2Jacobian*,void**); 456 457 PETSC_EXTERN PetscErrorCode TSRHSSplitSetIS(TS,const char[],IS); 458 PETSC_EXTERN PetscErrorCode TSRHSSplitGetIS(TS,const char[],IS*); 459 PETSC_EXTERN PetscErrorCode TSRHSSplitSetRHSFunction(TS,const char[],Vec,TSRHSFunction,void*); 460 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTS(TS,const char[],TS*); 461 PETSC_EXTERN PetscErrorCode TSRHSSplitGetSubTSs(TS,PetscInt*,TS*[]); 462 PETSC_EXTERN PetscErrorCode TSSetUseSplitRHSFunction(TS, PetscBool); 463 PETSC_EXTERN PetscErrorCode TSGetUseSplitRHSFunction(TS, PetscBool*); 464 465 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 466 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*); 467 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 468 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 469 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 470 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec); 471 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 472 473 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 474 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 475 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*)); 476 PETSC_EXTERN PetscErrorCode TSSetPostEvaluate(TS, PetscErrorCode(*)(TS)); 477 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 478 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 479 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 480 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*); 481 PETSC_EXTERN PetscErrorCode TSPostEvaluate(TS); 482 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 483 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 484 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 485 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 486 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*); 487 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*,PetscReal*,PetscReal*); 488 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*); 489 PETSC_EXTERN PetscErrorCode TSErrorWeightedENormInfinity(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*); 490 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm2(TS,Vec,Vec,Vec,PetscReal*,PetscReal*,PetscReal*); 491 PETSC_EXTERN PetscErrorCode TSErrorWeightedENorm(TS,Vec,Vec,Vec,NormType,PetscReal*,PetscReal*,PetscReal*); 492 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 493 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 494 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*)); 495 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*); 496 497 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 498 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*); 499 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 500 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 501 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 502 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *); 503 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 504 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 505 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 506 507 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 508 509 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 510 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat); 511 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 512 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool); 513 PETSC_EXTERN PetscErrorCode TSComputeI2Function(TS,PetscReal,Vec,Vec,Vec,Vec); 514 PETSC_EXTERN PetscErrorCode TSComputeI2Jacobian(TS,PetscReal,Vec,Vec,Vec,PetscReal,PetscReal,Mat,Mat); 515 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 516 517 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 518 519 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 520 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 521 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 522 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 523 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 524 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 525 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 526 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 527 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 528 PETSC_EXTERN PetscErrorCode DMTSSetI2Function(DM,TSI2Function,void*); 529 PETSC_EXTERN PetscErrorCode DMTSGetI2Function(DM,TSI2Function*,void**); 530 PETSC_EXTERN PetscErrorCode DMTSSetI2Jacobian(DM,TSI2Jacobian,void*); 531 PETSC_EXTERN PetscErrorCode DMTSGetI2Jacobian(DM,TSI2Jacobian*,void**); 532 533 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 534 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 535 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,TSForcingFunction,void*); 536 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,TSForcingFunction*,void**); 537 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*); 538 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal); 539 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **); 540 541 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*); 542 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*); 543 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*); 544 545 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 546 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 547 548 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 549 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*); 550 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 551 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*); 552 553 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 554 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *); 555 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 556 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *); 557 558 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM,Vec*,Vec*,PetscReal*); 559 PETSC_EXTERN PetscErrorCode DMPlexTSGetGradientDM(DM,PetscFV,DM*); 560 561 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 562 typedef struct { 563 Vec ray; 564 VecScatter scatter; 565 PetscViewer viewer; 566 TSMonitorLGCtx lgctx; 567 } TSMonitorDMDARayCtx; 568 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**); 569 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*); 570 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*); 571 572 /* Dynamic creation and loading functions */ 573 PETSC_EXTERN PetscFunctionList TSList; 574 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 575 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 576 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 577 578 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 579 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 580 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 581 582 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 583 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 584 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 585 PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 586 587 #define TS_FILE_CLASSID 1211225 588 589 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 590 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 591 592 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 593 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 594 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 595 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 596 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*); 597 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **); 598 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *); 599 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*); 600 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*); 601 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*); 602 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *); 603 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 604 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 605 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 606 PETSC_EXTERN PetscErrorCode TSMonitorError(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat *); 607 608 typedef struct _n_TSMonitorEnvelopeCtx* TSMonitorEnvelopeCtx; 609 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*); 610 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*); 611 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*); 612 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*); 613 614 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 615 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 616 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 617 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 618 619 typedef struct _n_TSMonitorSPCtx* TSMonitorSPCtx; 620 PETSC_EXTERN PetscErrorCode TSMonitorSPCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPCtx*); 621 PETSC_EXTERN PetscErrorCode TSMonitorSPCtxDestroy(TSMonitorSPCtx*); 622 PETSC_EXTERN PetscErrorCode TSMonitorSPSwarmSolution(TS,PetscInt,PetscReal,Vec,void*); 623 624 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*); 625 PETSC_EXTERN PetscErrorCode TSSetPostEventIntervalStep(TS,PetscReal); 626 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]); 627 /*J 628 TSSSPType - string with the name of TSSSP scheme. 629 630 Level: beginner 631 632 .seealso: TSSSPSetType(), TS 633 J*/ 634 typedef const char* TSSSPType; 635 #define TSSSPRKS2 "rks2" 636 #define TSSSPRKS3 "rks3" 637 #define TSSSPRK104 "rk104" 638 639 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 640 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 641 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 642 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 643 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 644 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 645 PETSC_EXTERN PetscFunctionList TSSSPList; 646 647 /*S 648 TSAdapt - Abstract object that manages time-step adaptivity 649 650 Level: beginner 651 652 .seealso: TS, TSAdaptCreate(), TSAdaptType 653 S*/ 654 typedef struct _p_TSAdapt *TSAdapt; 655 656 /*E 657 TSAdaptType - String with the name of TSAdapt scheme. 658 659 Level: beginner 660 661 .seealso: TSAdaptSetType(), TS 662 E*/ 663 typedef const char *TSAdaptType; 664 #define TSADAPTNONE "none" 665 #define TSADAPTBASIC "basic" 666 #define TSADAPTDSP "dsp" 667 #define TSADAPTCFL "cfl" 668 #define TSADAPTGLEE "glee" 669 #define TSADAPTHISTORY "history" 670 671 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*); 672 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt)); 673 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 674 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 675 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 676 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 677 PETSC_EXTERN PetscErrorCode TSAdaptGetType(TSAdapt,TSAdaptType*); 678 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 679 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 680 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 681 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 682 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 683 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*); 684 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 685 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 686 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt); 687 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt); 688 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 689 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 690 PETSC_EXTERN PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt,PetscBool); 691 PETSC_EXTERN PetscErrorCode TSAdaptSetSafety(TSAdapt,PetscReal,PetscReal); 692 PETSC_EXTERN PetscErrorCode TSAdaptGetSafety(TSAdapt,PetscReal*,PetscReal*); 693 PETSC_EXTERN PetscErrorCode TSAdaptSetClip(TSAdapt,PetscReal,PetscReal); 694 PETSC_EXTERN PetscErrorCode TSAdaptGetClip(TSAdapt,PetscReal*,PetscReal*); 695 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 696 PETSC_EXTERN PetscErrorCode TSAdaptGetStepLimits(TSAdapt,PetscReal*,PetscReal*); 697 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*)); 698 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetHistory(TSAdapt,PetscInt n,PetscReal hist[],PetscBool); 699 PETSC_EXTERN PetscErrorCode TSAdaptHistorySetTrajectory(TSAdapt,TSTrajectory,PetscBool); 700 PETSC_EXTERN PetscErrorCode TSAdaptHistoryGetStep(TSAdapt,PetscInt,PetscReal*,PetscReal*); 701 PETSC_EXTERN PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt,PetscInt); 702 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetFilter(TSAdapt,const char *); 703 PETSC_EXTERN PetscErrorCode TSAdaptDSPSetPID(TSAdapt,PetscReal,PetscReal,PetscReal); 704 705 /*S 706 TSGLLEAdapt - Abstract object that manages time-step adaptivity 707 708 Level: beginner 709 710 Developer Notes: 711 This functionality should be replaced by the TSAdapt. 712 713 .seealso: TSGLLE, TSGLLEAdaptCreate(), TSGLLEAdaptType 714 S*/ 715 typedef struct _p_TSGLLEAdapt *TSGLLEAdapt; 716 717 /*J 718 TSGLLEAdaptType - String with the name of TSGLLEAdapt scheme 719 720 Level: beginner 721 722 .seealso: TSGLLEAdaptSetType(), TS 723 J*/ 724 typedef const char *TSGLLEAdaptType; 725 #define TSGLLEADAPT_NONE "none" 726 #define TSGLLEADAPT_SIZE "size" 727 #define TSGLLEADAPT_BOTH "both" 728 729 PETSC_EXTERN PetscErrorCode TSGLLEAdaptRegister(const char[],PetscErrorCode (*)(TSGLLEAdapt)); 730 PETSC_EXTERN PetscErrorCode TSGLLEAdaptInitializePackage(void); 731 PETSC_EXTERN PetscErrorCode TSGLLEAdaptFinalizePackage(void); 732 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate(MPI_Comm,TSGLLEAdapt*); 733 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetType(TSGLLEAdapt,TSGLLEAdaptType); 734 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt,const char[]); 735 PETSC_EXTERN PetscErrorCode TSGLLEAdaptChoose(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 736 PETSC_EXTERN PetscErrorCode TSGLLEAdaptView(TSGLLEAdapt,PetscViewer); 737 PETSC_EXTERN PetscErrorCode TSGLLEAdaptSetFromOptions(PetscOptionItems*,TSGLLEAdapt); 738 PETSC_EXTERN PetscErrorCode TSGLLEAdaptDestroy(TSGLLEAdapt*); 739 740 /*J 741 TSGLLEAcceptType - String with the name of TSGLLEAccept scheme 742 743 Level: beginner 744 745 .seealso: TSGLLESetAcceptType(), TS 746 J*/ 747 typedef const char *TSGLLEAcceptType; 748 #define TSGLLEACCEPT_ALWAYS "always" 749 750 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLLEAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 751 PETSC_EXTERN PetscErrorCode TSGLLEAcceptRegister(const char[],TSGLLEAcceptFunction); 752 753 /*J 754 TSGLLEType - family of time integration method within the General Linear class 755 756 Level: beginner 757 758 .seealso: TSGLLESetType(), TSGLLERegister() 759 J*/ 760 typedef const char* TSGLLEType; 761 #define TSGLLE_IRKS "irks" 762 763 PETSC_EXTERN PetscErrorCode TSGLLERegister(const char[],PetscErrorCode(*)(TS)); 764 PETSC_EXTERN PetscErrorCode TSGLLEInitializePackage(void); 765 PETSC_EXTERN PetscErrorCode TSGLLEFinalizePackage(void); 766 PETSC_EXTERN PetscErrorCode TSGLLESetType(TS,TSGLLEType); 767 PETSC_EXTERN PetscErrorCode TSGLLEGetAdapt(TS,TSGLLEAdapt*); 768 PETSC_EXTERN PetscErrorCode TSGLLESetAcceptType(TS,TSGLLEAcceptType); 769 770 /*J 771 TSEIMEXType - String with the name of an Extrapolated IMEX method. 772 773 Level: beginner 774 775 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister() 776 J*/ 777 #define TSEIMEXType char* 778 779 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt); 780 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt); 781 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool); 782 783 /*J 784 TSRKType - String with the name of a Runge-Kutta method. 785 786 Level: beginner 787 788 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister() 789 J*/ 790 typedef const char* TSRKType; 791 #define TSRK1FE "1fe" 792 #define TSRK2A "2a" 793 #define TSRK3 "3" 794 #define TSRK3BS "3bs" 795 #define TSRK4 "4" 796 #define TSRK5F "5f" 797 #define TSRK5DP "5dp" 798 #define TSRK5BS "5bs" 799 800 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*); 801 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType); 802 PETSC_EXTERN PetscErrorCode TSRKSetMultirate(TS,PetscBool); 803 PETSC_EXTERN PetscErrorCode TSRKGetMultirate(TS,PetscBool*); 804 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool); 805 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 806 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 807 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 808 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 809 810 /*J 811 TSMPRKType - String with the name of a Partitioned Runge-Kutta method 812 813 Level: beginner 814 815 .seealso: TSMPRKSetType(), TS, TSMPRK, TSMPRKRegister() 816 J*/ 817 typedef const char* TSMPRKType; 818 #define TSMPRK2A22 "2a22" 819 #define TSMPRK2A23 "2a23" 820 #define TSMPRK2A32 "2a32" 821 #define TSMPRK2A33 "2a33" 822 #define TSMPRKP2 "p2" 823 #define TSMPRKP3 "p3" 824 825 PETSC_EXTERN PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType*); 826 PETSC_EXTERN PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType); 827 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[]); 828 PETSC_EXTERN PetscErrorCode TSMPRKInitializePackage(void); 829 PETSC_EXTERN PetscErrorCode TSMPRKFinalizePackage(void); 830 PETSC_EXTERN PetscErrorCode TSMPRKRegisterDestroy(void); 831 832 /*J 833 TSGLEEType - String with the name of a General Linear with Error Estimation method. 834 835 Level: beginner 836 837 .seealso: TSGLEESetType(), TS, TSGLEE, TSGLEERegister() 838 J*/ 839 typedef const char* TSGLEEType; 840 #define TSGLEEi1 "BE1" 841 #define TSGLEE23 "23" 842 #define TSGLEE24 "24" 843 #define TSGLEE25I "25i" 844 #define TSGLEE35 "35" 845 #define TSGLEEEXRK2A "exrk2a" 846 #define TSGLEERK32G1 "rk32g1" 847 #define TSGLEERK285EX "rk285ex" 848 /*J 849 TSGLEEMode - String with the mode of error estimation for a General Linear with Error Estimation method. 850 851 Level: beginner 852 853 .seealso: TSGLEESetMode(), TS, TSGLEE, TSGLEERegister() 854 J*/ 855 PETSC_EXTERN PetscErrorCode TSGLEEGetType(TS ts,TSGLEEType*); 856 PETSC_EXTERN PetscErrorCode TSGLEESetType(TS ts,TSGLEEType); 857 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[]); 858 PETSC_EXTERN PetscErrorCode TSGLEEFinalizePackage(void); 859 PETSC_EXTERN PetscErrorCode TSGLEEInitializePackage(void); 860 PETSC_EXTERN PetscErrorCode TSGLEERegisterDestroy(void); 861 862 /*J 863 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 864 865 Level: beginner 866 867 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 868 J*/ 869 typedef const char* TSARKIMEXType; 870 #define TSARKIMEX1BEE "1bee" 871 #define TSARKIMEXA2 "a2" 872 #define TSARKIMEXL2 "l2" 873 #define TSARKIMEXARS122 "ars122" 874 #define TSARKIMEX2C "2c" 875 #define TSARKIMEX2D "2d" 876 #define TSARKIMEX2E "2e" 877 #define TSARKIMEXPRSSP2 "prssp2" 878 #define TSARKIMEX3 "3" 879 #define TSARKIMEXBPR3 "bpr3" 880 #define TSARKIMEXARS443 "ars443" 881 #define TSARKIMEX4 "4" 882 #define TSARKIMEX5 "5" 883 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 884 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 885 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 886 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[]); 887 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 888 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 889 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 890 891 /*J 892 TSRosWType - String with the name of a Rosenbrock-W method. 893 894 Level: beginner 895 896 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 897 J*/ 898 typedef const char* TSRosWType; 899 #define TSROSW2M "2m" 900 #define TSROSW2P "2p" 901 #define TSROSWRA3PW "ra3pw" 902 #define TSROSWRA34PW2 "ra34pw2" 903 #define TSROSWRODAS3 "rodas3" 904 #define TSROSWSANDU3 "sandu3" 905 #define TSROSWASSP3P3S1C "assp3p3s1c" 906 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 907 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 908 #define TSROSWARK3 "ark3" 909 #define TSROSWTHETA1 "theta1" 910 #define TSROSWTHETA2 "theta2" 911 #define TSROSWGRK4T "grk4t" 912 #define TSROSWSHAMP4 "shamp4" 913 #define TSROSWVELDD4 "veldd4" 914 #define TSROSW4L "4l" 915 916 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS,TSRosWType*); 917 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS,TSRosWType); 918 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 919 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 920 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 921 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 922 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 923 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 924 925 PETSC_EXTERN PetscErrorCode TSBDFSetOrder(TS,PetscInt); 926 PETSC_EXTERN PetscErrorCode TSBDFGetOrder(TS,PetscInt*); 927 928 /*J 929 TSBasicSymplecticType - String with the name of a basic symplectic integration method. 930 931 Level: beginner 932 933 .seealso: TSBasicSymplecticSetType(), TS, TSBASICSYMPLECTIC, TSBasicSymplecticRegister() 934 J*/ 935 typedef const char* TSBasicSymplecticType; 936 #define TSBASICSYMPLECTICSIEULER "1" 937 #define TSBASICSYMPLECTICVELVERLET "2" 938 #define TSBASICSYMPLECTIC3 "3" 939 #define TSBASICSYMPLECTIC4 "4" 940 PETSC_EXTERN PetscErrorCode TSBasicSymplecticSetType(TS,TSBasicSymplecticType); 941 PETSC_EXTERN PetscErrorCode TSBasicSymplecticGetType(TS,TSBasicSymplecticType*); 942 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegister(TSBasicSymplecticType,PetscInt,PetscInt,PetscReal[],PetscReal[]); 943 PETSC_EXTERN PetscErrorCode TSBasicSymplecticInitializePackage(void); 944 PETSC_EXTERN PetscErrorCode TSBasicSymplecticFinalizePackage(void); 945 PETSC_EXTERN PetscErrorCode TSBasicSymplecticRegisterDestroy(void); 946 947 /* 948 PETSc interface to Sundials 949 */ 950 #ifdef PETSC_HAVE_SUNDIALS 951 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 952 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 953 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 954 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 955 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 956 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 957 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 958 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 959 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 960 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 961 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 962 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 963 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 964 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 965 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 966 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 967 #endif 968 969 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 970 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 971 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 972 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 973 974 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 975 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 976 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 977 978 PETSC_EXTERN PetscErrorCode TSAlpha2SetRadius(TS,PetscReal); 979 PETSC_EXTERN PetscErrorCode TSAlpha2SetParams(TS,PetscReal,PetscReal,PetscReal,PetscReal); 980 PETSC_EXTERN PetscErrorCode TSAlpha2GetParams(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*); 981 982 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 983 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 984 985 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 986 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*); 987 988 PETSC_EXTERN PetscErrorCode TSRHSJacobianTest(TS,PetscBool*); 989 PETSC_EXTERN PetscErrorCode TSRHSJacobianTestTranspose(TS,PetscBool*); 990 #endif 991