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 TSPSEUDO "pseudo" 31 #define TSCN "cn" 32 #define TSSUNDIALS "sundials" 33 #define TSRK "rk" 34 #define TSPYTHON "python" 35 #define TSTHETA "theta" 36 #define TSALPHA "alpha" 37 #define TSGL "gl" 38 #define TSSSP "ssp" 39 #define TSARKIMEX "arkimex" 40 #define TSROSW "rosw" 41 #define TSEIMEX "eimex" 42 #define TSMIMEX "mimex" 43 /*E 44 TSProblemType - Determines the type of problem this TS object is to be used to solve 45 46 Level: beginner 47 48 .seealso: TSCreate() 49 E*/ 50 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType; 51 52 /*E 53 TSEquationType - type of TS problem that is solved 54 55 Level: beginner 56 57 Developer Notes: this must match petsc/finclude/petscts.h 58 59 Supported types are: 60 TS_EQ_UNSPECIFIED (default) 61 TS_EQ_EXPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) := M(t) U_t - G(U,t) = 0 62 TS_EQ_IMPLICIT {ODE and DAE index 1, 2, 3, HI}: F(t,U,U_t) = 0 63 64 .seealso: TSGetEquationType(), TSSetEquationType() 65 E*/ 66 typedef enum { 67 TS_EQ_UNSPECIFIED = -1, 68 TS_EQ_EXPLICIT = 0, 69 TS_EQ_ODE_EXPLICIT = 1, 70 TS_EQ_DAE_SEMI_EXPLICIT_INDEX1 = 100, 71 TS_EQ_DAE_SEMI_EXPLICIT_INDEX2 = 200, 72 TS_EQ_DAE_SEMI_EXPLICIT_INDEX3 = 300, 73 TS_EQ_DAE_SEMI_EXPLICIT_INDEXHI = 500, 74 TS_EQ_IMPLICIT = 1000, 75 TS_EQ_ODE_IMPLICIT = 1001, 76 TS_EQ_DAE_IMPLICIT_INDEX1 = 1100, 77 TS_EQ_DAE_IMPLICIT_INDEX2 = 1200, 78 TS_EQ_DAE_IMPLICIT_INDEX3 = 1300, 79 TS_EQ_DAE_IMPLICIT_INDEXHI = 1500 80 } TSEquationType; 81 PETSC_EXTERN const char *const*TSEquationTypes; 82 83 /*E 84 TSConvergedReason - reason a TS method has converged or not 85 86 Level: beginner 87 88 Developer Notes: this must match petsc/finclude/petscts.h 89 90 Each reason has its own manual page. 91 92 .seealso: TSGetConvergedReason() 93 E*/ 94 typedef enum { 95 TS_CONVERGED_ITERATING = 0, 96 TS_CONVERGED_TIME = 1, 97 TS_CONVERGED_ITS = 2, 98 TS_CONVERGED_USER = 3, 99 TS_CONVERGED_EVENT = 4, 100 TS_CONVERGED_PSEUDO_FATOL = 5, 101 TS_CONVERGED_PSEUDO_FRTOL = 6, 102 TS_DIVERGED_NONLINEAR_SOLVE = -1, 103 TS_DIVERGED_STEP_REJECTED = -2 104 } TSConvergedReason; 105 PETSC_EXTERN const char *const*TSConvergedReasons; 106 107 /*MC 108 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 109 110 Level: beginner 111 112 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt() 113 M*/ 114 115 /*MC 116 TS_CONVERGED_TIME - the final time was reached 117 118 Level: beginner 119 120 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime() 121 M*/ 122 123 /*MC 124 TS_CONVERGED_ITS - the maximum number of iterations (time-steps) was reached prior to the final time 125 126 Level: beginner 127 128 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetDuration() 129 M*/ 130 131 /*MC 132 TS_CONVERGED_USER - user requested termination 133 134 Level: beginner 135 136 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 137 M*/ 138 139 /*MC 140 TS_CONVERGED_EVENT - user requested termination on event detection 141 142 Level: beginner 143 144 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration() 145 M*/ 146 147 /*MC 148 TS_CONVERGED_PSEUDO_FRTOL - stops when function norm decreased by a set amount, used only for TSPSEUDO 149 150 Level: beginner 151 152 Options Database: 153 . -ts_pseudo_frtol <rtol> 154 155 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration(), TS_CONVERGED_PSEUDO_FATOL 156 M*/ 157 158 /*MC 159 TS_CONVERGED_PSEUDO_FATOL - stops when function norm decreases below a set amount, used only for TSPSEUDO 160 161 Level: beginner 162 163 Options Database: 164 . -ts_pseudo_fatol <atol> 165 166 .seealso: TSSolve(), TSGetConvergedReason(), TSSetConvergedReason(), TSSetDuration(), TS_CONVERGED_PSEUDO_FRTOL 167 M*/ 168 169 /*MC 170 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 171 172 Level: beginner 173 174 Notes: See TSSetMaxSNESFailures() for how to allow more nonlinear solver failures. 175 176 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason(), TSSetMaxSNESFailures() 177 M*/ 178 179 /*MC 180 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 181 182 Level: beginner 183 184 Notes: See TSSetMaxStepRejections() for how to allow more step rejections. 185 186 .seealso: TSSolve(), TSGetConvergedReason(), TSGetAdapt(), TSSetMaxStepRejections() 187 M*/ 188 189 /*E 190 TSExactFinalTimeOption - option for handling of final time step 191 192 Level: beginner 193 194 Developer Notes: this must match petsc/finclude/petscts.h 195 196 $ TS_EXACTFINALTIME_STEPOVER - Don't do anything if final time is exceeded 197 $ TS_EXACTFINALTIME_INTERPOLATE - Interpolate back to final time 198 $ TS_EXACTFINALTIME_MATCHSTEP - Adapt final time step to match the final time 199 200 .seealso: TSGetConvergedReason(), TSSetExactFinalTime() 201 202 E*/ 203 typedef enum {TS_EXACTFINALTIME_UNSPECIFIED=0,TS_EXACTFINALTIME_STEPOVER=1,TS_EXACTFINALTIME_INTERPOLATE=2,TS_EXACTFINALTIME_MATCHSTEP=3} TSExactFinalTimeOption; 204 PETSC_EXTERN const char *const TSExactFinalTimeOptions[]; 205 206 207 /* Logging support */ 208 PETSC_EXTERN PetscClassId TS_CLASSID; 209 PETSC_EXTERN PetscClassId DMTS_CLASSID; 210 PETSC_EXTERN PetscClassId TSADAPT_CLASSID; 211 212 PETSC_EXTERN PetscErrorCode TSInitializePackage(void); 213 214 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 215 PETSC_EXTERN PetscErrorCode TSClone(TS,TS*); 216 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 217 218 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 219 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 220 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 221 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 222 PETSC_EXTERN PetscErrorCode TSMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*)); 223 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 224 225 PETSC_EXTERN PetscErrorCode TSAdjointMonitor(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*); 226 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *,PetscErrorCode (*)(void**)); 227 PETSC_EXTERN PetscErrorCode TSAdjointMonitorCancel(TS); 228 PETSC_EXTERN PetscErrorCode TSAdjointMonitorSetFromOptions(TS,const char[],const char[],const char[],PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*)(TS,PetscViewerAndFormat*)); 229 230 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 231 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 232 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 233 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 234 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 235 PETSC_EXTERN PetscErrorCode TSReset(TS); 236 237 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 238 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 239 240 /*S 241 TSTrajectory - Abstract PETSc object that storing the trajectory (solution of ODE/ADE at each time step and stage) 242 243 Level: advanced 244 245 Concepts: ODE solvers, adjoint 246 247 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC, TSDestroy() 248 S*/ 249 typedef struct _p_TSTrajectory* TSTrajectory; 250 251 /*J 252 TSTrajectoryType - String with the name of a PETSc TS trajectory storage method 253 254 Level: intermediate 255 256 .seealso: TSSetType(), TS, TSRegister(), TSTrajectoryCreate(), TSTrajectorySetType() 257 J*/ 258 typedef const char* TSTrajectoryType; 259 #define TSTRAJECTORYBASIC "basic" 260 #define TSTRAJECTORYSINGLEFILE "singlefile" 261 #define TSTRAJECTORYMEMORY "memory" 262 #define TSTRAJECTORYVISUALIZATION "visualization" 263 264 PETSC_EXTERN PetscFunctionList TSTrajectoryList; 265 PETSC_EXTERN PetscClassId TSTRAJECTORY_CLASSID; 266 PETSC_EXTERN PetscBool TSTrajectoryRegisterAllCalled; 267 268 PETSC_EXTERN PetscErrorCode TSSetSaveTrajectory(TS); 269 270 PETSC_EXTERN PetscErrorCode TSTrajectoryCreate(MPI_Comm,TSTrajectory*); 271 PETSC_EXTERN PetscErrorCode TSTrajectoryDestroy(TSTrajectory*); 272 PETSC_EXTERN PetscErrorCode TSTrajectorySetType(TSTrajectory,TS,const TSTrajectoryType); 273 PETSC_EXTERN PetscErrorCode TSTrajectorySet(TSTrajectory,TS,PetscInt,PetscReal,Vec); 274 PETSC_EXTERN PetscErrorCode TSTrajectoryGet(TSTrajectory,TS,PetscInt,PetscReal*); 275 PETSC_EXTERN PetscErrorCode TSTrajectorySetFromOptions(TSTrajectory,TS); 276 PETSC_EXTERN PetscErrorCode TSTrajectoryRegisterAll(void); 277 PETSC_EXTERN PetscErrorCode TSTrajectorySetUp(TSTrajectory,TS); 278 279 PETSC_EXTERN PetscErrorCode TSSetCostGradients(TS,PetscInt,Vec*,Vec*); 280 PETSC_EXTERN PetscErrorCode TSGetCostGradients(TS,PetscInt*,Vec**,Vec**); 281 PETSC_EXTERN PetscErrorCode TSSetCostIntegrand(TS,PetscInt,PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*),PetscBool,void*); 282 PETSC_EXTERN PetscErrorCode TSGetCostIntegral(TS,Vec*); 283 284 PETSC_EXTERN PetscErrorCode TSAdjointSetRHSJacobian(TS,Mat,PetscErrorCode(*)(TS,PetscReal,Vec,Mat,void*),void*); 285 PETSC_EXTERN PetscErrorCode TSAdjointSolve(TS); 286 PETSC_EXTERN PetscErrorCode TSAdjointSetSteps(TS,PetscInt); 287 288 PETSC_EXTERN PetscErrorCode TSAdjointComputeRHSJacobian(TS,PetscReal,Vec,Mat); 289 PETSC_EXTERN PetscErrorCode TSAdjointStep(TS); 290 PETSC_EXTERN PetscErrorCode TSAdjointSetUp(TS); 291 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDPFunction(TS,PetscReal,Vec,Vec*); 292 PETSC_EXTERN PetscErrorCode TSAdjointComputeDRDYFunction(TS,PetscReal,Vec,Vec*); 293 PETSC_EXTERN PetscErrorCode TSAdjointComputeCostIntegrand(TS,PetscReal,Vec,Vec); 294 PETSC_EXTERN PetscErrorCode TSAdjointCostIntegral(TS); 295 PETSC_EXTERN PetscErrorCode TSForwardCostIntegral(TS); 296 297 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 298 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 299 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,TSExactFinalTimeOption); 300 301 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*); 302 303 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 304 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 305 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 306 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 307 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolutionPhase(TS,PetscInt,PetscReal,Vec,void*); 308 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 309 310 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDefault(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat *); 311 PETSC_EXTERN PetscErrorCode TSAdjointMonitorDrawSensi(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*); 312 313 PETSC_EXTERN PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,PetscViewerAndFormat*); 314 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 315 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 316 317 PETSC_EXTERN PetscErrorCode TSStep(TS); 318 PETSC_EXTERN PetscErrorCode TSEvaluateWLTE(TS,NormType,PetscInt*,PetscReal*); 319 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 320 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 321 PETSC_EXTERN PetscErrorCode TSGetEquationType(TS,TSEquationType*); 322 PETSC_EXTERN PetscErrorCode TSSetEquationType(TS,TSEquationType); 323 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 324 PETSC_EXTERN PetscErrorCode TSSetConvergedReason(TS,TSConvergedReason); 325 PETSC_EXTERN PetscErrorCode TSGetSolveTime(TS,PetscReal*); 326 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 327 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 328 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 329 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 330 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 331 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 332 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 333 PETSC_EXTERN PetscErrorCode TSRollBack(TS); 334 PETSC_EXTERN PetscErrorCode TSGetTotalSteps(TS,PetscInt*); 335 336 PETSC_EXTERN PetscErrorCode TSGetStages(TS,PetscInt*,Vec**); 337 338 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 339 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 340 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 341 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 342 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 343 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 344 PETSC_EXTERN PetscErrorCode TSGetPrevTime(TS,PetscReal*); 345 346 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 347 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat,Mat,void*); 348 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 349 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 350 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 351 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 352 PETSC_EXTERN PetscErrorCode TSRHSJacobianSetReuse(TS,PetscBool); 353 354 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 355 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 356 PETSC_EXTERN PetscErrorCode TSSetForcingFunction(TS,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*); 357 358 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 359 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 360 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 361 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 362 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 363 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 364 365 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 366 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat,Mat,void*); 367 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 368 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 369 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 370 PETSC_EXTERN PetscErrorCode TSComputeForcingFunction(TS,PetscReal,Vec); 371 PETSC_EXTERN PetscErrorCode TSComputeIJacobianDefaultColor(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 372 373 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 374 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 375 PETSC_EXTERN PetscErrorCode TSSetPostStage(TS, PetscErrorCode (*)(TS,PetscReal,PetscInt,Vec*)); 376 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 377 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 378 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 379 PETSC_EXTERN PetscErrorCode TSPostStage(TS,PetscReal,PetscInt,Vec*); 380 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 381 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 382 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 383 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 384 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 385 PETSC_EXTERN PetscErrorCode TSErrorWeightedNormInfinity(TS,Vec,Vec,PetscReal*); 386 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm2(TS,Vec,Vec,PetscReal*); 387 PETSC_EXTERN PetscErrorCode TSErrorWeightedNorm(TS,Vec,Vec,NormType,PetscReal*); 388 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 389 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 390 PETSC_EXTERN PetscErrorCode TSSetFunctionDomainError(TS, PetscErrorCode (*)(TS,PetscReal,Vec,PetscBool*)); 391 PETSC_EXTERN PetscErrorCode TSFunctionDomainError(TS,PetscReal,Vec,PetscBool*); 392 393 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 394 PETSC_EXTERN PetscErrorCode TSPseudoTimeStepDefault(TS,PetscReal*,void*); 395 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 396 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 397 398 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 399 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStepDefault(TS,Vec,void*,PetscReal*,PetscBool *); 400 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 401 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 402 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 403 404 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 405 406 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 407 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat,Mat); 408 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 409 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,PetscBool); 410 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 411 412 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 413 414 PETSC_EXTERN PetscErrorCode DMTSSetBoundaryLocal(DM, PetscErrorCode (*)(DM, PetscReal, Vec, Vec, void *), void *); 415 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 416 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 417 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 418 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 419 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 420 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 421 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 422 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 423 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 424 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 425 PETSC_EXTERN PetscErrorCode DMTSSetForcingFunction(DM,PetscErrorCode (*)(TS,PetscReal,Vec,void*),void*); 426 PETSC_EXTERN PetscErrorCode DMTSGetForcingFunction(DM,PetscErrorCode (**)(TS,PetscReal,Vec,void*),void**); 427 PETSC_EXTERN PetscErrorCode DMTSGetMinRadius(DM,PetscReal*); 428 PETSC_EXTERN PetscErrorCode DMTSSetMinRadius(DM,PetscReal); 429 PETSC_EXTERN PetscErrorCode DMTSCheckFromOptions(TS, Vec, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void **); 430 431 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,Vec,void*),void*); 432 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*),void*); 433 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunctionLocal(DM,PetscErrorCode (*)(DM,PetscReal,Vec,Vec,void*),void*); 434 435 PETSC_EXTERN PetscErrorCode DMTSSetIFunctionSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 436 PETSC_EXTERN PetscErrorCode DMTSSetIJacobianSerialize(DM,PetscErrorCode (*)(void*,PetscViewer),PetscErrorCode (*)(void**,PetscViewer)); 437 438 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 439 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*); 440 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 441 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*); 442 443 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 444 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,void*),void *); 445 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 446 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*),void *); 447 448 PETSC_EXTERN PetscErrorCode DMPlexTSGetGeometryFVM(DM, Vec *, Vec *, PetscReal *); 449 450 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 451 typedef struct { 452 Vec ray; 453 VecScatter scatter; 454 PetscViewer viewer; 455 TSMonitorLGCtx lgctx; 456 } TSMonitorDMDARayCtx; 457 PETSC_EXTERN PetscErrorCode TSMonitorDMDARayDestroy(void**); 458 PETSC_EXTERN PetscErrorCode TSMonitorDMDARay(TS,PetscInt,PetscReal,Vec,void*); 459 PETSC_EXTERN PetscErrorCode TSMonitorLGDMDARay(TS,PetscInt,PetscReal,Vec,void*); 460 461 462 /* Dynamic creation and loading functions */ 463 PETSC_EXTERN PetscFunctionList TSList; 464 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 465 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 466 PETSC_EXTERN PetscErrorCode TSRegister(const char[], PetscErrorCode (*)(TS)); 467 468 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 469 PETSC_EXTERN PetscErrorCode TSSetSNES(TS,SNES); 470 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 471 472 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 473 PETSC_EXTERN PetscErrorCode TSLoad(TS,PetscViewer); 474 PETSC_STATIC_INLINE PetscErrorCode TSViewFromOptions(TS A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 475 PETSC_STATIC_INLINE PetscErrorCode TSTrajectoryViewFromOptions(TSTrajectory A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 476 477 #define TS_FILE_CLASSID 1211225 478 479 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 480 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 481 482 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 483 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 484 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 485 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 486 PETSC_EXTERN PetscErrorCode TSMonitorLGSetVariableNames(TS,const char * const*); 487 PETSC_EXTERN PetscErrorCode TSMonitorLGGetVariableNames(TS,const char *const **); 488 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetVariableNames(TSMonitorLGCtx,const char * const *); 489 PETSC_EXTERN PetscErrorCode TSMonitorLGSetDisplayVariables(TS,const char * const*); 490 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetDisplayVariables(TSMonitorLGCtx,const char * const*); 491 PETSC_EXTERN PetscErrorCode TSMonitorLGSetTransform(TS,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void*); 492 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxSetTransform(TSMonitorLGCtx,PetscErrorCode (*)(void*,Vec,Vec*),PetscErrorCode (*)(void*),void *); 493 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 494 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 495 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 496 497 typedef struct _n_TSMonitorEnvelopeCtx* TSMonitorEnvelopeCtx; 498 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxCreate(TS,TSMonitorEnvelopeCtx*); 499 PETSC_EXTERN PetscErrorCode TSMonitorEnvelope(TS,PetscInt,PetscReal,Vec,void*); 500 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeGetBounds(TS,Vec*,Vec*); 501 PETSC_EXTERN PetscErrorCode TSMonitorEnvelopeCtxDestroy(TSMonitorEnvelopeCtx*); 502 503 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 504 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 505 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 506 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 507 508 PETSC_EXTERN PetscErrorCode TSSetEventHandler(TS,PetscInt,PetscInt[],PetscBool[],PetscErrorCode (*)(TS,PetscReal,Vec,PetscScalar[],void*),PetscErrorCode (*)(TS,PetscInt,PetscInt[],PetscReal,Vec,PetscBool,void*),void*); 509 PETSC_EXTERN PetscErrorCode TSSetEventTolerances(TS,PetscReal,PetscReal[]); 510 /*J 511 TSSSPType - string with the name of TSSSP scheme. 512 513 Level: beginner 514 515 .seealso: TSSSPSetType(), TS 516 J*/ 517 typedef const char* TSSSPType; 518 #define TSSSPRKS2 "rks2" 519 #define TSSSPRKS3 "rks3" 520 #define TSSSPRK104 "rk104" 521 522 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 523 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 524 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 525 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 526 PETSC_EXTERN PetscErrorCode TSSSPFinalizePackage(void); 527 PETSC_EXTERN PetscErrorCode TSSSPInitializePackage(void); 528 529 /*S 530 TSAdapt - Abstract object that manages time-step adaptivity 531 532 Level: beginner 533 534 .seealso: TS, TSAdaptCreate(), TSAdaptType 535 S*/ 536 typedef struct _p_TSAdapt *TSAdapt; 537 538 /*E 539 TSAdaptType - String with the name of TSAdapt scheme. 540 541 Level: beginner 542 543 .seealso: TSAdaptSetType(), TS 544 E*/ 545 typedef const char *TSAdaptType; 546 #define TSADAPTBASIC "basic" 547 #define TSADAPTNONE "none" 548 #define TSADAPTCFL "cfl" 549 550 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*); 551 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],PetscErrorCode (*)(TSAdapt)); 552 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(void); 553 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 554 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 555 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 556 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 557 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 558 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 559 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 560 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 561 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscReal,Vec,PetscBool*); 562 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 563 PETSC_EXTERN PetscErrorCode TSAdaptLoad(TSAdapt,PetscViewer); 564 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems*,TSAdapt); 565 PETSC_EXTERN PetscErrorCode TSAdaptReset(TSAdapt); 566 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 567 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 568 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 569 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscReal,Vec,PetscBool*)); 570 571 PETSC_EXTERN PetscErrorCode TSAdaptBasicSetClip(TSAdapt,PetscReal,PetscReal); 572 PETSC_EXTERN PetscErrorCode TSAdaptBasicGetClip(TSAdapt,PetscReal*,PetscReal*); 573 574 /*S 575 TSGLAdapt - Abstract object that manages time-step adaptivity 576 577 Level: beginner 578 579 Developer Notes: 580 This functionality should be replaced by the TSAdapt. 581 582 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 583 S*/ 584 typedef struct _p_TSGLAdapt *TSGLAdapt; 585 586 /*J 587 TSGLAdaptType - String with the name of TSGLAdapt scheme 588 589 Level: beginner 590 591 .seealso: TSGLAdaptSetType(), TS 592 J*/ 593 typedef const char *TSGLAdaptType; 594 #define TSGLADAPT_NONE "none" 595 #define TSGLADAPT_SIZE "size" 596 #define TSGLADAPT_BOTH "both" 597 598 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],PetscErrorCode (*)(TSGLAdapt)); 599 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(void); 600 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 601 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 602 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 603 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 604 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 605 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 606 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(PetscOptionItems*,TSGLAdapt); 607 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 608 609 /*J 610 TSGLAcceptType - String with the name of TSGLAccept scheme 611 612 Level: beginner 613 614 .seealso: TSGLSetAcceptType(), TS 615 J*/ 616 typedef const char *TSGLAcceptType; 617 #define TSGLACCEPT_ALWAYS "always" 618 619 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 620 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],TSGLAcceptFunction); 621 622 /*J 623 TSGLType - family of time integration method within the General Linear class 624 625 Level: beginner 626 627 .seealso: TSGLSetType(), TSGLRegister() 628 J*/ 629 typedef const char* TSGLType; 630 #define TSGL_IRKS "irks" 631 632 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],PetscErrorCode(*)(TS)); 633 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(void); 634 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 635 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 636 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 637 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 638 639 /*J 640 TSEIMEXType - String with the name of an Extrapolated IMEX method. 641 642 Level: beginner 643 644 .seealso: TSEIMEXSetType(), TS, TSEIMEX, TSEIMEXRegister() 645 J*/ 646 #define TSEIMEXType char* 647 648 PETSC_EXTERN PetscErrorCode TSEIMEXSetMaxRows(TS ts,PetscInt); 649 PETSC_EXTERN PetscErrorCode TSEIMEXSetRowCol(TS ts,PetscInt,PetscInt); 650 PETSC_EXTERN PetscErrorCode TSEIMEXSetOrdAdapt(TS,PetscBool); 651 652 /*J 653 TSRKType - String with the name of a Runge-Kutta method. 654 655 Level: beginner 656 657 .seealso: TSRKSetType(), TS, TSRK, TSRKRegister() 658 J*/ 659 typedef const char* TSRKType; 660 #define TSRK1FE "1fe" 661 #define TSRK2A "2a" 662 #define TSRK3 "3" 663 #define TSRK3BS "3bs" 664 #define TSRK4 "4" 665 #define TSRK5F "5f" 666 #define TSRK5DP "5dp" 667 PETSC_EXTERN PetscErrorCode TSRKGetType(TS ts,TSRKType*); 668 PETSC_EXTERN PetscErrorCode TSRKSetType(TS ts,TSRKType); 669 PETSC_EXTERN PetscErrorCode TSRKSetFullyImplicit(TS,PetscBool); 670 PETSC_EXTERN PetscErrorCode TSRKRegister(TSRKType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 671 PETSC_EXTERN PetscErrorCode TSRKFinalizePackage(void); 672 PETSC_EXTERN PetscErrorCode TSRKInitializePackage(void); 673 PETSC_EXTERN PetscErrorCode TSRKRegisterDestroy(void); 674 675 /*J 676 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 677 678 Level: beginner 679 680 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 681 J*/ 682 typedef const char* TSARKIMEXType; 683 #define TSARKIMEX1BEE "1bee" 684 #define TSARKIMEXA2 "a2" 685 #define TSARKIMEXL2 "l2" 686 #define TSARKIMEXARS122 "ars122" 687 #define TSARKIMEX2C "2c" 688 #define TSARKIMEX2D "2d" 689 #define TSARKIMEX2E "2e" 690 #define TSARKIMEXPRSSP2 "prssp2" 691 #define TSARKIMEX3 "3" 692 #define TSARKIMEXBPR3 "bpr3" 693 #define TSARKIMEXARS443 "ars443" 694 #define TSARKIMEX4 "4" 695 #define TSARKIMEX5 "5" 696 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 697 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 698 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 699 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[]); 700 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 701 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(void); 702 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 703 704 /*J 705 TSRosWType - String with the name of a Rosenbrock-W method. 706 707 Level: beginner 708 709 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 710 J*/ 711 typedef const char* TSRosWType; 712 #define TSROSW2M "2m" 713 #define TSROSW2P "2p" 714 #define TSROSWRA3PW "ra3pw" 715 #define TSROSWRA34PW2 "ra34pw2" 716 #define TSROSWRODAS3 "rodas3" 717 #define TSROSWSANDU3 "sandu3" 718 #define TSROSWASSP3P3S1C "assp3p3s1c" 719 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 720 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 721 #define TSROSWARK3 "ark3" 722 #define TSROSWTHETA1 "theta1" 723 #define TSROSWTHETA2 "theta2" 724 #define TSROSWGRK4T "grk4t" 725 #define TSROSWSHAMP4 "shamp4" 726 #define TSROSWVELDD4 "veldd4" 727 #define TSROSW4L "4l" 728 729 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 730 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 731 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 732 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 733 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 734 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 735 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(void); 736 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 737 738 /* 739 PETSc interface to Sundials 740 */ 741 #ifdef PETSC_HAVE_SUNDIALS 742 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 743 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 744 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 745 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 746 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 747 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 748 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 749 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 750 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 751 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 752 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 753 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 754 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 755 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 756 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 757 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 758 #endif 759 760 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 761 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 762 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 763 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 764 765 PETSC_EXTERN PetscErrorCode TSAlphaUseAdapt(TS,PetscBool); 766 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 767 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 768 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 769 770 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 771 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 772 773 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 774 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat,Mat,void*); 775 776 #endif 777