1 /* 2 User interface for the timestepping package. This package 3 is for use in solving time-dependent PDEs. 4 */ 5 #if !defined(__PETSCTS_H) 6 #define __PETSCTS_H 7 #include <petscsnes.h> 8 9 /*S 10 TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 11 12 Level: beginner 13 14 Concepts: ODE solvers 15 16 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC 17 S*/ 18 typedef struct _p_TS* TS; 19 20 /*J 21 TSType - String with the name of a PETSc TS method or the creation function 22 with an optional dynamic library name, for example 23 http://www.mcs.anl.gov/petsc/lib.a:mytscreate() 24 25 Level: beginner 26 27 .seealso: TSSetType(), TS 28 J*/ 29 typedef const char* TSType; 30 #define TSEULER "euler" 31 #define TSBEULER "beuler" 32 #define TSPSEUDO "pseudo" 33 #define TSCN "cn" 34 #define TSSUNDIALS "sundials" 35 #define TSRK "rk" 36 #define TSPYTHON "python" 37 #define TSTHETA "theta" 38 #define TSALPHA "alpha" 39 #define TSGL "gl" 40 #define TSSSP "ssp" 41 #define TSARKIMEX "arkimex" 42 #define TSROSW "rosw" 43 44 /*E 45 TSProblemType - Determines the type of problem this TS object is to be used to solve 46 47 Level: beginner 48 49 .seealso: TSCreate() 50 E*/ 51 typedef enum {TS_LINEAR,TS_NONLINEAR} TSProblemType; 52 53 /*E 54 TSConvergedReason - reason a TS method has converged or not 55 56 Level: beginner 57 58 Developer Notes: this must match finclude/petscts.h 59 60 Each reason has its own manual page. 61 62 .seealso: TSGetConvergedReason() 63 E*/ 64 typedef enum { 65 TS_CONVERGED_ITERATING = 0, 66 TS_CONVERGED_TIME = 1, 67 TS_CONVERGED_ITS = 2, 68 TS_DIVERGED_NONLINEAR_SOLVE = -1, 69 TS_DIVERGED_STEP_REJECTED = -2 70 } TSConvergedReason; 71 PETSC_EXTERN const char *const*TSConvergedReasons; 72 73 /*MC 74 TS_CONVERGED_ITERATING - this only occurs if TSGetConvergedReason() is called during the TSSolve() 75 76 Level: beginner 77 78 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt() 79 M*/ 80 81 /*MC 82 TS_CONVERGED_TIME - the final time was reached 83 84 Level: beginner 85 86 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration(), TSGetSolveTime() 87 M*/ 88 89 /*MC 90 TS_CONVERGED_ITS - the maximum number of iterations was reached prior to the final time 91 92 Level: beginner 93 94 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSSetDuration() 95 M*/ 96 97 /*MC 98 TS_DIVERGED_NONLINEAR_SOLVE - too many nonlinear solves failed 99 100 Level: beginner 101 102 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt(), TSGetSNES(), SNESGetConvergedReason() 103 M*/ 104 105 /*MC 106 TS_DIVERGED_STEP_REJECTED - too many steps were rejected 107 108 Level: beginner 109 110 .seealso: TSSolve(), TSConvergedReason(), TSGetAdapt() 111 M*/ 112 113 /* Logging support */ 114 PETSC_EXTERN PetscClassId TS_CLASSID; 115 116 PETSC_EXTERN PetscErrorCode TSInitializePackage(const char[]); 117 118 PETSC_EXTERN PetscErrorCode TSCreate(MPI_Comm,TS*); 119 PETSC_EXTERN PetscErrorCode TSDestroy(TS*); 120 121 PETSC_EXTERN PetscErrorCode TSSetProblemType(TS,TSProblemType); 122 PETSC_EXTERN PetscErrorCode TSGetProblemType(TS,TSProblemType*); 123 PETSC_EXTERN PetscErrorCode TSMonitor(TS,PetscInt,PetscReal,Vec); 124 PETSC_EXTERN PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 125 PETSC_EXTERN PetscErrorCode TSMonitorCancel(TS); 126 127 PETSC_EXTERN PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 128 PETSC_EXTERN PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 129 PETSC_EXTERN PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 130 PETSC_EXTERN PetscErrorCode TSSetFromOptions(TS); 131 PETSC_EXTERN PetscErrorCode TSSetUp(TS); 132 PETSC_EXTERN PetscErrorCode TSReset(TS); 133 134 PETSC_EXTERN PetscErrorCode TSSetSolution(TS,Vec); 135 PETSC_EXTERN PetscErrorCode TSGetSolution(TS,Vec*); 136 137 PETSC_EXTERN PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 138 PETSC_EXTERN PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 139 PETSC_EXTERN PetscErrorCode TSSetExactFinalTime(TS,PetscBool); 140 141 PETSC_EXTERN PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 142 143 typedef struct _n_TSMonitorDrawCtx* TSMonitorDrawCtx; 144 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorDrawCtx *); 145 PETSC_EXTERN PetscErrorCode TSMonitorDrawCtxDestroy(TSMonitorDrawCtx*); 146 PETSC_EXTERN PetscErrorCode TSMonitorDrawSolution(TS,PetscInt,PetscReal,Vec,void*); 147 PETSC_EXTERN PetscErrorCode TSMonitorDrawError(TS,PetscInt,PetscReal,Vec,void*); 148 149 PETSC_EXTERN PetscErrorCode TSMonitorSolutionBinary(TS,PetscInt,PetscReal,Vec,void*); 150 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTK(TS,PetscInt,PetscReal,Vec,void*); 151 PETSC_EXTERN PetscErrorCode TSMonitorSolutionVTKDestroy(void*); 152 153 PETSC_EXTERN PetscErrorCode TSStep(TS); 154 PETSC_EXTERN PetscErrorCode TSEvaluateStep(TS,PetscInt,Vec,PetscBool*); 155 PETSC_EXTERN PetscErrorCode TSSolve(TS,Vec); 156 PETSC_EXTERN PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 157 PETSC_EXTERN PetscErrorCode TSGetSolutionTime(TS,PetscReal*); 158 PETSC_EXTERN PetscErrorCode TSGetSNESIterations(TS,PetscInt*); 159 PETSC_EXTERN PetscErrorCode TSGetKSPIterations(TS,PetscInt*); 160 PETSC_EXTERN PetscErrorCode TSGetStepRejections(TS,PetscInt*); 161 PETSC_EXTERN PetscErrorCode TSSetMaxStepRejections(TS,PetscInt); 162 PETSC_EXTERN PetscErrorCode TSGetSNESFailures(TS,PetscInt*); 163 PETSC_EXTERN PetscErrorCode TSSetMaxSNESFailures(TS,PetscInt); 164 PETSC_EXTERN PetscErrorCode TSSetErrorIfStepFails(TS,PetscBool); 165 166 PETSC_EXTERN PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 167 PETSC_EXTERN PetscErrorCode TSGetTimeStep(TS,PetscReal*); 168 PETSC_EXTERN PetscErrorCode TSGetTime(TS,PetscReal*); 169 PETSC_EXTERN PetscErrorCode TSSetTime(TS,PetscReal); 170 PETSC_EXTERN PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 171 PETSC_EXTERN PetscErrorCode TSSetTimeStep(TS,PetscReal); 172 173 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 174 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 175 PETSC_EXTERN PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 176 PETSC_EXTERN PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 177 PETSC_EXTERN PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 178 PETSC_EXTERN PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 179 180 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSSolutionFunction)(TS,PetscReal,Vec,void*); 181 PETSC_EXTERN PetscErrorCode TSSetSolutionFunction(TS,TSSolutionFunction,void*); 182 183 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 184 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 185 PETSC_EXTERN PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 186 PETSC_EXTERN PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 187 PETSC_EXTERN PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 188 PETSC_EXTERN PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 189 190 PETSC_EXTERN PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 191 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 192 PETSC_EXTERN PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 193 PETSC_EXTERN PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 194 PETSC_EXTERN PetscErrorCode TSComputeSolutionFunction(TS,PetscReal,Vec); 195 196 PETSC_EXTERN PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 197 PETSC_EXTERN PetscErrorCode TSSetPreStage(TS, PetscErrorCode (*)(TS,PetscReal)); 198 PETSC_EXTERN PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 199 PETSC_EXTERN PetscErrorCode TSPreStep(TS); 200 PETSC_EXTERN PetscErrorCode TSPreStage(TS,PetscReal); 201 PETSC_EXTERN PetscErrorCode TSPostStep(TS); 202 PETSC_EXTERN PetscErrorCode TSSetRetainStages(TS,PetscBool); 203 PETSC_EXTERN PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 204 PETSC_EXTERN PetscErrorCode TSSetTolerances(TS,PetscReal,Vec,PetscReal,Vec); 205 PETSC_EXTERN PetscErrorCode TSGetTolerances(TS,PetscReal*,Vec*,PetscReal*,Vec*); 206 PETSC_EXTERN PetscErrorCode TSErrorNormWRMS(TS,Vec,PetscReal*); 207 PETSC_EXTERN PetscErrorCode TSSetCFLTimeLocal(TS,PetscReal); 208 PETSC_EXTERN PetscErrorCode TSGetCFLTime(TS,PetscReal*); 209 210 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 211 PETSC_EXTERN PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 212 PETSC_EXTERN PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 213 PETSC_EXTERN PetscErrorCode TSPseudoSetMaxTimeStep(TS,PetscReal); 214 215 PETSC_EXTERN PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 216 PETSC_EXTERN PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 217 PETSC_EXTERN PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 218 PETSC_EXTERN PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 219 PETSC_EXTERN PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 220 221 PETSC_EXTERN PetscErrorCode TSPythonSetType(TS,const char[]); 222 223 PETSC_EXTERN PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 224 PETSC_EXTERN PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 225 PETSC_EXTERN PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 226 PETSC_EXTERN PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 227 PETSC_EXTERN PetscErrorCode TSComputeLinearStability(TS,PetscReal,PetscReal,PetscReal*,PetscReal*); 228 229 PETSC_EXTERN PetscErrorCode TSVISetVariableBounds(TS,Vec,Vec); 230 231 PETSC_EXTERN PetscErrorCode DMTSSetRHSFunction(DM,TSRHSFunction,void*); 232 PETSC_EXTERN PetscErrorCode DMTSGetRHSFunction(DM,TSRHSFunction*,void**); 233 PETSC_EXTERN PetscErrorCode DMTSSetRHSJacobian(DM,TSRHSJacobian,void*); 234 PETSC_EXTERN PetscErrorCode DMTSGetRHSJacobian(DM,TSRHSJacobian*,void**); 235 PETSC_EXTERN PetscErrorCode DMTSSetIFunction(DM,TSIFunction,void*); 236 PETSC_EXTERN PetscErrorCode DMTSGetIFunction(DM,TSIFunction*,void**); 237 PETSC_EXTERN PetscErrorCode DMTSSetIJacobian(DM,TSIJacobian,void*); 238 PETSC_EXTERN PetscErrorCode DMTSGetIJacobian(DM,TSIJacobian*,void**); 239 PETSC_EXTERN PetscErrorCode DMTSSetSolutionFunction(DM,TSSolutionFunction,void*); 240 PETSC_EXTERN PetscErrorCode DMTSGetSolutionFunction(DM,TSSolutionFunction*,void**); 241 242 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 243 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSRHSJacobianLocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*); 244 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIFunctionLocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 245 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*DMDATSIJacobianLocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*); 246 247 PETSC_EXTERN PetscErrorCode DMDATSSetRHSFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*),void *); 248 PETSC_EXTERN PetscErrorCode DMDATSSetRHSJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*),void *); 249 PETSC_EXTERN PetscErrorCode DMDATSSetIFunctionLocal(DM,InsertMode,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*),void *); 250 PETSC_EXTERN PetscErrorCode DMDATSSetIJacobianLocal(DM,PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*),void *); 251 252 /* Dynamic creation and loading functions */ 253 PETSC_EXTERN PetscFList TSList; 254 PETSC_EXTERN PetscBool TSRegisterAllCalled; 255 PETSC_EXTERN PetscErrorCode TSGetType(TS,TSType*); 256 PETSC_EXTERN PetscErrorCode TSSetType(TS,TSType); 257 PETSC_EXTERN PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 258 PETSC_EXTERN PetscErrorCode TSRegisterAll(const char[]); 259 PETSC_EXTERN PetscErrorCode TSRegisterDestroy(void); 260 261 /*MC 262 TSRegisterDynamic - Adds a creation method to the TS package. 263 264 Synopsis: 265 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 266 267 Not Collective 268 269 Input Parameters: 270 + name - The name of a new user-defined creation routine 271 . path - The path (either absolute or relative) of the library containing this routine 272 . func_name - The name of the creation routine 273 - create_func - The creation routine itself 274 275 Notes: 276 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 277 278 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 279 280 Sample usage: 281 .vb 282 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 283 .ve 284 285 Then, your ts type can be chosen with the procedural interface via 286 .vb 287 TS ts; 288 TSCreate(MPI_Comm, &ts); 289 TSSetType(ts, "my_ts") 290 .ve 291 or at runtime via the option 292 .vb 293 -ts_type my_ts 294 .ve 295 296 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 297 If your function is not being put into a shared library then use TSRegister() instead 298 299 Level: advanced 300 301 .keywords: TS, register 302 .seealso: TSRegisterAll(), TSRegisterDestroy() 303 M*/ 304 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 305 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 306 #else 307 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 308 #endif 309 310 PETSC_EXTERN PetscErrorCode TSGetSNES(TS,SNES*); 311 PETSC_EXTERN PetscErrorCode TSGetKSP(TS,KSP*); 312 313 PETSC_EXTERN PetscErrorCode TSView(TS,PetscViewer); 314 315 PETSC_EXTERN PetscErrorCode TSSetApplicationContext(TS,void *); 316 PETSC_EXTERN PetscErrorCode TSGetApplicationContext(TS,void *); 317 318 typedef struct _n_TSMonitorLGCtx* TSMonitorLGCtx; 319 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorLGCtx *); 320 PETSC_EXTERN PetscErrorCode TSMonitorLGCtxDestroy(TSMonitorLGCtx*); 321 PETSC_EXTERN PetscErrorCode TSMonitorLGTimeStep(TS,PetscInt,PetscReal,Vec,void *); 322 PETSC_EXTERN PetscErrorCode TSMonitorLGSolution(TS,PetscInt,PetscReal,Vec,void *); 323 PETSC_EXTERN PetscErrorCode TSMonitorLGError(TS,PetscInt,PetscReal,Vec,void *); 324 PETSC_EXTERN PetscErrorCode TSMonitorLGSNESIterations(TS,PetscInt,PetscReal,Vec,void *); 325 PETSC_EXTERN PetscErrorCode TSMonitorLGKSPIterations(TS,PetscInt,PetscReal,Vec,void *); 326 327 typedef struct _n_TSMonitorSPEigCtx* TSMonitorSPEigCtx; 328 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscInt,TSMonitorSPEigCtx *); 329 PETSC_EXTERN PetscErrorCode TSMonitorSPEigCtxDestroy(TSMonitorSPEigCtx*); 330 PETSC_EXTERN PetscErrorCode TSMonitorSPEig(TS,PetscInt,PetscReal,Vec,void *); 331 332 /*J 333 TSSSPType - string with the name of TSSSP scheme. 334 335 Level: beginner 336 337 .seealso: TSSSPSetType(), TS 338 J*/ 339 typedef const char* TSSSPType; 340 #define TSSSPRKS2 "rks2" 341 #define TSSSPRKS3 "rks3" 342 #define TSSSPRK104 "rk104" 343 344 PETSC_EXTERN PetscErrorCode TSSSPSetType(TS,TSSSPType); 345 PETSC_EXTERN PetscErrorCode TSSSPGetType(TS,TSSSPType*); 346 PETSC_EXTERN PetscErrorCode TSSSPSetNumStages(TS,PetscInt); 347 PETSC_EXTERN PetscErrorCode TSSSPGetNumStages(TS,PetscInt*); 348 349 /*S 350 TSAdapt - Abstract object that manages time-step adaptivity 351 352 Level: beginner 353 354 .seealso: TS, TSAdaptCreate(), TSAdaptType 355 S*/ 356 typedef struct _p_TSAdapt *TSAdapt; 357 358 /*E 359 TSAdaptType - String with the name of TSAdapt scheme or the creation function 360 with an optional dynamic library name, for example 361 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 362 363 Level: beginner 364 365 .seealso: TSAdaptSetType(), TS 366 E*/ 367 typedef const char *TSAdaptType; 368 #define TSADAPTBASIC "basic" 369 #define TSADAPTNONE "none" 370 #define TSADAPTCFL "cfl" 371 372 /*MC 373 TSAdaptRegisterDynamic - adds a TSAdapt implementation 374 375 Synopsis: 376 PetscErrorCode TSAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 377 378 Not Collective 379 380 Input Parameters: 381 + name_scheme - name of user-defined adaptivity scheme 382 . path - path (either absolute or relative) the library containing this scheme 383 . name_create - name of routine to create method context 384 - routine_create - routine to create method context 385 386 Notes: 387 TSAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 388 389 If dynamic libraries are used, then the fourth input argument (routine_create) 390 is ignored. 391 392 Sample usage: 393 .vb 394 TSAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 395 "MySchemeCreate",MySchemeCreate); 396 .ve 397 398 Then, your scheme can be chosen with the procedural interface via 399 $ TSAdaptSetType(ts,"my_scheme") 400 or at runtime via the option 401 $ -ts_adapt_type my_scheme 402 403 Level: advanced 404 405 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 406 and others of the form ${any_environmental_variable} occuring in pathname will be 407 replaced with appropriate values. 408 409 .keywords: TSAdapt, register 410 411 .seealso: TSAdaptRegisterAll() 412 M*/ 413 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 414 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,0) 415 #else 416 # define TSAdaptRegisterDynamic(a,b,c,d) TSAdaptRegister(a,b,c,d) 417 #endif 418 419 PETSC_EXTERN PetscErrorCode TSGetAdapt(TS,TSAdapt*); 420 PETSC_EXTERN PetscErrorCode TSAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSAdapt)); 421 PETSC_EXTERN PetscErrorCode TSAdaptRegisterAll(const char[]); 422 PETSC_EXTERN PetscErrorCode TSAdaptRegisterDestroy(void); 423 PETSC_EXTERN PetscErrorCode TSAdaptInitializePackage(const char[]); 424 PETSC_EXTERN PetscErrorCode TSAdaptFinalizePackage(void); 425 PETSC_EXTERN PetscErrorCode TSAdaptCreate(MPI_Comm,TSAdapt*); 426 PETSC_EXTERN PetscErrorCode TSAdaptSetType(TSAdapt,TSAdaptType); 427 PETSC_EXTERN PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt,const char[]); 428 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesClear(TSAdapt); 429 PETSC_EXTERN PetscErrorCode TSAdaptCandidateAdd(TSAdapt,const char[],PetscInt,PetscInt,PetscReal,PetscReal,PetscBool); 430 PETSC_EXTERN PetscErrorCode TSAdaptCandidatesGet(TSAdapt,PetscInt*,const PetscInt**,const PetscInt**,const PetscReal**,const PetscReal**); 431 PETSC_EXTERN PetscErrorCode TSAdaptChoose(TSAdapt,TS,PetscReal,PetscInt*,PetscReal*,PetscBool*); 432 PETSC_EXTERN PetscErrorCode TSAdaptCheckStage(TSAdapt,TS,PetscBool*); 433 PETSC_EXTERN PetscErrorCode TSAdaptView(TSAdapt,PetscViewer); 434 PETSC_EXTERN PetscErrorCode TSAdaptSetFromOptions(TSAdapt); 435 PETSC_EXTERN PetscErrorCode TSAdaptDestroy(TSAdapt*); 436 PETSC_EXTERN PetscErrorCode TSAdaptSetMonitor(TSAdapt,PetscBool); 437 PETSC_EXTERN PetscErrorCode TSAdaptSetStepLimits(TSAdapt,PetscReal,PetscReal); 438 PETSC_EXTERN PetscErrorCode TSAdaptSetCheckStage(TSAdapt,PetscErrorCode(*)(TSAdapt,TS,PetscBool*)); 439 440 /*S 441 TSGLAdapt - Abstract object that manages time-step adaptivity 442 443 Level: beginner 444 445 Developer Notes: 446 This functionality should be replaced by the TSAdapt. 447 448 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 449 S*/ 450 typedef struct _p_TSGLAdapt *TSGLAdapt; 451 452 /*J 453 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 454 with an optional dynamic library name, for example 455 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 456 457 Level: beginner 458 459 .seealso: TSGLAdaptSetType(), TS 460 J*/ 461 typedef const char *TSGLAdaptType; 462 #define TSGLADAPT_NONE "none" 463 #define TSGLADAPT_SIZE "size" 464 #define TSGLADAPT_BOTH "both" 465 466 /*MC 467 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 468 469 Synopsis: 470 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 471 472 Not Collective 473 474 Input Parameters: 475 + name_scheme - name of user-defined adaptivity scheme 476 . path - path (either absolute or relative) the library containing this scheme 477 . name_create - name of routine to create method context 478 - routine_create - routine to create method context 479 480 Notes: 481 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 482 483 If dynamic libraries are used, then the fourth input argument (routine_create) 484 is ignored. 485 486 Sample usage: 487 .vb 488 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 489 "MySchemeCreate",MySchemeCreate); 490 .ve 491 492 Then, your scheme can be chosen with the procedural interface via 493 $ TSGLAdaptSetType(ts,"my_scheme") 494 or at runtime via the option 495 $ -ts_adapt_type my_scheme 496 497 Level: advanced 498 499 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 500 and others of the form ${any_environmental_variable} occuring in pathname will be 501 replaced with appropriate values. 502 503 .keywords: TSGLAdapt, register 504 505 .seealso: TSGLAdaptRegisterAll() 506 M*/ 507 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 508 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 509 #else 510 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 511 #endif 512 513 PETSC_EXTERN PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 514 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterAll(const char[]); 515 PETSC_EXTERN PetscErrorCode TSGLAdaptRegisterDestroy(void); 516 PETSC_EXTERN PetscErrorCode TSGLAdaptInitializePackage(const char[]); 517 PETSC_EXTERN PetscErrorCode TSGLAdaptFinalizePackage(void); 518 PETSC_EXTERN PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 519 PETSC_EXTERN PetscErrorCode TSGLAdaptSetType(TSGLAdapt,TSGLAdaptType); 520 PETSC_EXTERN PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 521 PETSC_EXTERN PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 522 PETSC_EXTERN PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 523 PETSC_EXTERN PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 524 PETSC_EXTERN PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 525 526 /*J 527 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 528 with an optional dynamic library name, for example 529 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 530 531 Level: beginner 532 533 .seealso: TSGLSetAcceptType(), TS 534 J*/ 535 typedef const char *TSGLAcceptType; 536 #define TSGLACCEPT_ALWAYS "always" 537 538 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 539 PETSC_EXTERN PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 540 541 /*MC 542 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 543 544 Synopsis: 545 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 546 547 Not Collective 548 549 Input Parameters: 550 + name_scheme - name of user-defined acceptance scheme 551 . path - path (either absolute or relative) the library containing this scheme 552 . name_create - name of routine to create method context 553 - routine_create - routine to create method context 554 555 Notes: 556 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 557 558 If dynamic libraries are used, then the fourth input argument (routine_create) 559 is ignored. 560 561 Sample usage: 562 .vb 563 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 564 "MySchemeCreate",MySchemeCreate); 565 .ve 566 567 Then, your scheme can be chosen with the procedural interface via 568 $ TSGLSetAcceptType(ts,"my_scheme") 569 or at runtime via the option 570 $ -ts_gl_accept_type my_scheme 571 572 Level: advanced 573 574 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 575 and others of the form ${any_environmental_variable} occuring in pathname will be 576 replaced with appropriate values. 577 578 .keywords: TSGL, TSGLAcceptType, register 579 580 .seealso: TSGLRegisterAll() 581 M*/ 582 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 583 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 584 #else 585 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 586 #endif 587 588 /*J 589 TSGLType - family of time integration method within the General Linear class 590 591 Level: beginner 592 593 .seealso: TSGLSetType(), TSGLRegister() 594 J*/ 595 typedef const char* TSGLType; 596 #define TSGL_IRKS "irks" 597 598 /*MC 599 TSGLRegisterDynamic - adds a TSGL implementation 600 601 Synopsis: 602 PetscErrorCode TSGLRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 603 604 Not Collective 605 606 Input Parameters: 607 + name_scheme - name of user-defined general linear scheme 608 . path - path (either absolute or relative) the library containing this scheme 609 . name_create - name of routine to create method context 610 - routine_create - routine to create method context 611 612 Notes: 613 TSGLRegisterDynamic() may be called multiple times to add several user-defined families. 614 615 If dynamic libraries are used, then the fourth input argument (routine_create) 616 is ignored. 617 618 Sample usage: 619 .vb 620 TSGLRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 621 "MySchemeCreate",MySchemeCreate); 622 .ve 623 624 Then, your scheme can be chosen with the procedural interface via 625 $ TSGLSetType(ts,"my_scheme") 626 or at runtime via the option 627 $ -ts_gl_type my_scheme 628 629 Level: advanced 630 631 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 632 and others of the form ${any_environmental_variable} occuring in pathname will be 633 replaced with appropriate values. 634 635 .keywords: TSGL, register 636 637 .seealso: TSGLRegisterAll() 638 M*/ 639 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 640 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 641 #else 642 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 643 #endif 644 645 PETSC_EXTERN PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 646 PETSC_EXTERN PetscErrorCode TSGLRegisterAll(const char[]); 647 PETSC_EXTERN PetscErrorCode TSGLRegisterDestroy(void); 648 PETSC_EXTERN PetscErrorCode TSGLInitializePackage(const char[]); 649 PETSC_EXTERN PetscErrorCode TSGLFinalizePackage(void); 650 PETSC_EXTERN PetscErrorCode TSGLSetType(TS,TSGLType); 651 PETSC_EXTERN PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 652 PETSC_EXTERN PetscErrorCode TSGLSetAcceptType(TS,TSGLAcceptType); 653 654 /*J 655 TSARKIMEXType - String with the name of an Additive Runge-Kutta IMEX method. 656 657 Level: beginner 658 659 .seealso: TSARKIMEXSetType(), TS, TSARKIMEX, TSARKIMEXRegister() 660 J*/ 661 typedef const char* TSARKIMEXType; 662 #define TSARKIMEXA2 "a2" 663 #define TSARKIMEXL2 "l2" 664 #define TSARKIMEXARS122 "ars122" 665 #define TSARKIMEX2C "2c" 666 #define TSARKIMEX2D "2d" 667 #define TSARKIMEX2E "2e" 668 #define TSARKIMEXPRSSP2 "prssp2" 669 #define TSARKIMEX3 "3" 670 #define TSARKIMEXBPR3 "bpr3" 671 #define TSARKIMEXARS443 "ars443" 672 #define TSARKIMEX4 "4" 673 #define TSARKIMEX5 "5" 674 PETSC_EXTERN PetscErrorCode TSARKIMEXGetType(TS ts,TSARKIMEXType*); 675 PETSC_EXTERN PetscErrorCode TSARKIMEXSetType(TS ts,TSARKIMEXType); 676 PETSC_EXTERN PetscErrorCode TSARKIMEXSetFullyImplicit(TS,PetscBool); 677 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[]); 678 PETSC_EXTERN PetscErrorCode TSARKIMEXFinalizePackage(void); 679 PETSC_EXTERN PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 680 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterDestroy(void); 681 PETSC_EXTERN PetscErrorCode TSARKIMEXRegisterAll(void); 682 683 /*J 684 TSRosWType - String with the name of a Rosenbrock-W method. 685 686 Level: beginner 687 688 .seealso: TSRosWSetType(), TS, TSROSW, TSRosWRegister() 689 J*/ 690 typedef const char* TSRosWType; 691 #define TSROSW2M "2m" 692 #define TSROSW2P "2p" 693 #define TSROSWRA3PW "ra3pw" 694 #define TSROSWRA34PW2 "ra34pw2" 695 #define TSROSWRODAS3 "rodas3" 696 #define TSROSWSANDU3 "sandu3" 697 #define TSROSWASSP3P3S1C "assp3p3s1c" 698 #define TSROSWLASSP3P4S2C "lassp3p4s2c" 699 #define TSROSWLLSSP3P4S2C "llssp3p4s2c" 700 #define TSROSWARK3 "ark3" 701 #define TSROSWTHETA1 "theta1" 702 #define TSROSWTHETA2 "theta2" 703 #define TSROSWGRK4T "grk4t" 704 #define TSROSWSHAMP4 "shamp4" 705 #define TSROSWVELDD4 "veldd4" 706 #define TSROSW4L "4l" 707 708 PETSC_EXTERN PetscErrorCode TSRosWGetType(TS ts,TSRosWType*); 709 PETSC_EXTERN PetscErrorCode TSRosWSetType(TS ts,TSRosWType); 710 PETSC_EXTERN PetscErrorCode TSRosWSetRecomputeJacobian(TS,PetscBool); 711 PETSC_EXTERN PetscErrorCode TSRosWRegister(TSRosWType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[]); 712 PETSC_EXTERN PetscErrorCode TSRosWRegisterRos4(TSRosWType,PetscReal,PetscReal,PetscReal,PetscReal,PetscReal); 713 PETSC_EXTERN PetscErrorCode TSRosWFinalizePackage(void); 714 PETSC_EXTERN PetscErrorCode TSRosWInitializePackage(const char path[]); 715 PETSC_EXTERN PetscErrorCode TSRosWRegisterDestroy(void); 716 PETSC_EXTERN PetscErrorCode TSRosWRegisterAll(void); 717 718 /* 719 PETSc interface to Sundials 720 */ 721 #ifdef PETSC_HAVE_SUNDIALS 722 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 723 PETSC_EXTERN const char *const TSSundialsLmmTypes[]; 724 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 725 PETSC_EXTERN const char *const TSSundialsGramSchmidtTypes[]; 726 PETSC_EXTERN PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 727 PETSC_EXTERN PetscErrorCode TSSundialsGetPC(TS,PC*); 728 PETSC_EXTERN PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 729 PETSC_EXTERN PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 730 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 731 PETSC_EXTERN PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 732 PETSC_EXTERN PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 733 PETSC_EXTERN PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 734 PETSC_EXTERN PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 735 PETSC_EXTERN PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 736 PETSC_EXTERN PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 737 PETSC_EXTERN PetscErrorCode TSSundialsSetMaxl(TS,PetscInt); 738 #endif 739 740 PETSC_EXTERN PetscErrorCode TSRKSetTolerance(TS,PetscReal); 741 742 PETSC_EXTERN PetscErrorCode TSThetaSetTheta(TS,PetscReal); 743 PETSC_EXTERN PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 744 PETSC_EXTERN PetscErrorCode TSThetaGetEndpoint(TS,PetscBool*); 745 PETSC_EXTERN PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 746 747 PETSC_EXTERN PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 748 PETSC_EXTERN PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 749 PETSC_EXTERN PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 750 PETSC_EXTERN PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 751 PETSC_EXTERN PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 752 753 PETSC_EXTERN PetscErrorCode TSSetDM(TS,DM); 754 PETSC_EXTERN PetscErrorCode TSGetDM(TS,DM*); 755 756 PETSC_EXTERN PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 757 PETSC_EXTERN PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 758 759 #endif 760