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