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