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