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