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