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