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 PETSC_EXTERN_CXX_BEGIN 9 10 /*S 11 TS - Abstract PETSc object that manages all time-steppers (ODE integrators) 12 13 Level: beginner 14 15 Concepts: ODE solvers 16 17 .seealso: TSCreate(), TSSetType(), TSType, SNES, KSP, PC 18 S*/ 19 typedef struct _p_TS* TS; 20 21 /*E 22 TSType - String with the name of a PETSc TS method or the creation function 23 with an optional dynamic library name, for example 24 http://www.mcs.anl.gov/petsc/lib.a:mytscreate() 25 26 Level: beginner 27 28 .seealso: TSSetType(), TS 29 E*/ 30 #define TSType char* 31 #define TSEULER "euler" 32 #define TSBEULER "beuler" 33 #define TSPSEUDO "pseudo" 34 #define TSCN "cn" 35 #define TSSUNDIALS "sundials" 36 #define TSRK "rk" 37 #define TSPYTHON "python" 38 #define TSTHETA "theta" 39 #define TSALPHA "alpha" 40 #define TSGL "gl" 41 #define TSSSP "ssp" 42 #define TSARKIMEX "arkimex" 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 typedef enum { 54 TS_CONVERGED_ITERATING = 0, 55 TS_CONVERGED_TIME = 1, 56 TS_CONVERGED_ITS = 2, 57 TS_DIVERGED_NONLINEAR_SOLVE = -1, 58 TS_DIVERGED_STEP_REJECTED = -2 59 } TSConvergedReason; 60 extern const char *const*TSConvergedReasons; 61 62 /* Logging support */ 63 extern PetscClassId TS_CLASSID; 64 65 extern PetscErrorCode TSInitializePackage(const char[]); 66 67 extern PetscErrorCode TSCreate(MPI_Comm,TS*); 68 extern PetscErrorCode TSDestroy(TS*); 69 70 extern PetscErrorCode TSSetProblemType(TS,TSProblemType); 71 extern PetscErrorCode TSGetProblemType(TS,TSProblemType*); 72 extern PetscErrorCode TSMonitorSet(TS,PetscErrorCode(*)(TS,PetscInt,PetscReal,Vec,void*),void *,PetscErrorCode (*)(void**)); 73 extern PetscErrorCode TSMonitorCancel(TS); 74 75 extern PetscErrorCode TSSetOptionsPrefix(TS,const char[]); 76 extern PetscErrorCode TSAppendOptionsPrefix(TS,const char[]); 77 extern PetscErrorCode TSGetOptionsPrefix(TS,const char *[]); 78 extern PetscErrorCode TSSetFromOptions(TS); 79 extern PetscErrorCode TSSetUp(TS); 80 extern PetscErrorCode TSReset(TS); 81 82 extern PetscErrorCode TSSetSolution(TS,Vec); 83 extern PetscErrorCode TSGetSolution(TS,Vec*); 84 85 extern PetscErrorCode TSSetDuration(TS,PetscInt,PetscReal); 86 extern PetscErrorCode TSGetDuration(TS,PetscInt*,PetscReal*); 87 extern PetscErrorCode TSSetExactFinalTime(TS,PetscBool); 88 89 extern PetscErrorCode TSMonitorDefault(TS,PetscInt,PetscReal,Vec,void*); 90 extern PetscErrorCode TSMonitorSolution(TS,PetscInt,PetscReal,Vec,void*); 91 extern PetscErrorCode TSStep(TS); 92 extern PetscErrorCode TSSolve(TS,Vec,PetscReal*); 93 extern PetscErrorCode TSGetConvergedReason(TS,TSConvergedReason*); 94 95 extern PetscErrorCode TSSetInitialTimeStep(TS,PetscReal,PetscReal); 96 extern PetscErrorCode TSGetTimeStep(TS,PetscReal*); 97 extern PetscErrorCode TSGetTime(TS,PetscReal*); 98 extern PetscErrorCode TSSetTime(TS,PetscReal); 99 extern PetscErrorCode TSGetTimeStepNumber(TS,PetscInt*); 100 extern PetscErrorCode TSSetTimeStep(TS,PetscReal); 101 102 typedef PetscErrorCode (*TSRHSFunction)(TS,PetscReal,Vec,Vec,void*); 103 typedef PetscErrorCode (*TSRHSJacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 104 extern PetscErrorCode TSSetRHSFunction(TS,Vec,TSRHSFunction,void*); 105 extern PetscErrorCode TSGetRHSFunction(TS,Vec*,TSRHSFunction*,void**); 106 extern PetscErrorCode TSSetRHSJacobian(TS,Mat,Mat,TSRHSJacobian,void*); 107 extern PetscErrorCode TSGetRHSJacobian(TS,Mat*,Mat*,TSRHSJacobian*,void**); 108 109 typedef PetscErrorCode (*TSIFunction)(TS,PetscReal,Vec,Vec,Vec,void*); 110 typedef PetscErrorCode (*TSIJacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 111 extern PetscErrorCode TSSetIFunction(TS,Vec,TSIFunction,void*); 112 extern PetscErrorCode TSGetIFunction(TS,Vec*,TSIFunction*,void**); 113 extern PetscErrorCode TSSetIJacobian(TS,Mat,Mat,TSIJacobian,void*); 114 extern PetscErrorCode TSGetIJacobian(TS,Mat*,Mat*,TSIJacobian*,void**); 115 116 extern PetscErrorCode TSComputeRHSFunctionLinear(TS,PetscReal,Vec,Vec,void*); 117 extern PetscErrorCode TSComputeRHSJacobianConstant(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 118 extern PetscErrorCode TSComputeIFunctionLinear(TS,PetscReal,Vec,Vec,Vec,void*); 119 extern PetscErrorCode TSComputeIJacobianConstant(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*); 120 extern PetscErrorCode TSDefaultComputeJacobianColor(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 121 extern PetscErrorCode TSDefaultComputeJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*); 122 123 extern PetscErrorCode TSSetPreStep(TS, PetscErrorCode (*)(TS)); 124 extern PetscErrorCode TSSetPostStep(TS, PetscErrorCode (*)(TS)); 125 extern PetscErrorCode TSPreStep(TS); 126 extern PetscErrorCode TSPostStep(TS); 127 extern PetscErrorCode TSSetRetainStages(TS,PetscBool); 128 extern PetscErrorCode TSInterpolate(TS,PetscReal,Vec); 129 130 extern PetscErrorCode TSPseudoSetTimeStep(TS,PetscErrorCode(*)(TS,PetscReal*,void*),void*); 131 extern PetscErrorCode TSPseudoDefaultTimeStep(TS,PetscReal*,void*); 132 extern PetscErrorCode TSPseudoComputeTimeStep(TS,PetscReal *); 133 134 extern PetscErrorCode TSPseudoSetVerifyTimeStep(TS,PetscErrorCode(*)(TS,Vec,void*,PetscReal*,PetscBool *),void*); 135 extern PetscErrorCode TSPseudoDefaultVerifyTimeStep(TS,Vec,void*,PetscReal*,PetscBool *); 136 extern PetscErrorCode TSPseudoVerifyTimeStep(TS,Vec,PetscReal*,PetscBool *); 137 extern PetscErrorCode TSPseudoSetTimeStepIncrement(TS,PetscReal); 138 extern PetscErrorCode TSPseudoIncrementDtFromInitialDt(TS); 139 140 extern PetscErrorCode TSPythonSetType(TS,const char[]); 141 142 extern PetscErrorCode TSComputeRHSFunction(TS,PetscReal,Vec,Vec); 143 extern PetscErrorCode TSComputeRHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*); 144 extern PetscErrorCode TSComputeIFunction(TS,PetscReal,Vec,Vec,Vec,PetscBool); 145 extern PetscErrorCode TSComputeIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,PetscBool); 146 147 /* Dynamic creation and loading functions */ 148 extern PetscFList TSList; 149 extern PetscBool TSRegisterAllCalled; 150 extern PetscErrorCode TSGetType(TS,const TSType*); 151 extern PetscErrorCode TSSetType(TS,const TSType); 152 extern PetscErrorCode TSRegister(const char[], const char[], const char[], PetscErrorCode (*)(TS)); 153 extern PetscErrorCode TSRegisterAll(const char[]); 154 extern PetscErrorCode TSRegisterDestroy(void); 155 156 /*MC 157 TSRegisterDynamic - Adds a creation method to the TS package. 158 159 Synopsis: 160 PetscErrorCode TSRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(TS)) 161 162 Not Collective 163 164 Input Parameters: 165 + name - The name of a new user-defined creation routine 166 . path - The path (either absolute or relative) of the library containing this routine 167 . func_name - The name of the creation routine 168 - create_func - The creation routine itself 169 170 Notes: 171 TSRegisterDynamic() may be called multiple times to add several user-defined tses. 172 173 If dynamic libraries are used, then the fourth input argument (create_func) is ignored. 174 175 Sample usage: 176 .vb 177 TSRegisterDynamic("my_ts", "/home/username/my_lib/lib/libO/solaris/libmy.a", "MyTSCreate", MyTSCreate); 178 .ve 179 180 Then, your ts type can be chosen with the procedural interface via 181 .vb 182 TS ts; 183 TSCreate(MPI_Comm, &ts); 184 TSSetType(ts, "my_ts") 185 .ve 186 or at runtime via the option 187 .vb 188 -ts_type my_ts 189 .ve 190 191 Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 192 If your function is not being put into a shared library then use TSRegister() instead 193 194 Level: advanced 195 196 .keywords: TS, register 197 .seealso: TSRegisterAll(), TSRegisterDestroy() 198 M*/ 199 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 200 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,0) 201 #else 202 #define TSRegisterDynamic(a,b,c,d) TSRegister(a,b,c,d) 203 #endif 204 205 extern PetscErrorCode TSGetSNES(TS,SNES*); 206 extern PetscErrorCode TSGetKSP(TS,KSP*); 207 208 extern PetscErrorCode TSView(TS,PetscViewer); 209 210 extern PetscErrorCode TSSetApplicationContext(TS,void *); 211 extern PetscErrorCode TSGetApplicationContext(TS,void *); 212 213 extern PetscErrorCode TSMonitorLGCreate(const char[],const char[],int,int,int,int,PetscDrawLG *); 214 extern PetscErrorCode TSMonitorLG(TS,PetscInt,PetscReal,Vec,void *); 215 extern PetscErrorCode TSMonitorLGDestroy(PetscDrawLG*); 216 217 /*S 218 TSGLAdapt - Abstract object that manages time-step adaptivity 219 220 Level: beginner 221 222 .seealso: TSGL, TSGLAdaptCreate(), TSGLAdaptType 223 S*/ 224 typedef struct _p_TSGLAdapt *TSGLAdapt; 225 226 /*E 227 TSGLAdaptType - String with the name of TSGLAdapt scheme or the creation function 228 with an optional dynamic library name, for example 229 http://www.mcs.anl.gov/petsc/lib.a:mytsgladaptcreate() 230 231 Level: beginner 232 233 .seealso: TSGLAdaptSetType(), TS 234 E*/ 235 #define TSGLAdaptType char* 236 #define TSGLADAPT_NONE "none" 237 #define TSGLADAPT_SIZE "size" 238 #define TSGLADAPT_BOTH "both" 239 240 /*MC 241 TSGLAdaptRegisterDynamic - adds a TSGLAdapt implementation 242 243 Synopsis: 244 PetscErrorCode TSGLAdaptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 245 246 Not Collective 247 248 Input Parameters: 249 + name_scheme - name of user-defined adaptivity scheme 250 . path - path (either absolute or relative) the library containing this scheme 251 . name_create - name of routine to create method context 252 - routine_create - routine to create method context 253 254 Notes: 255 TSGLAdaptRegisterDynamic() may be called multiple times to add several user-defined families. 256 257 If dynamic libraries are used, then the fourth input argument (routine_create) 258 is ignored. 259 260 Sample usage: 261 .vb 262 TSGLAdaptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 263 "MySchemeCreate",MySchemeCreate); 264 .ve 265 266 Then, your scheme can be chosen with the procedural interface via 267 $ TSGLAdaptSetType(ts,"my_scheme") 268 or at runtime via the option 269 $ -ts_adapt_type my_scheme 270 271 Level: advanced 272 273 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 274 and others of the form ${any_environmental_variable} occuring in pathname will be 275 replaced with appropriate values. 276 277 .keywords: TSGLAdapt, register 278 279 .seealso: TSGLAdaptRegisterAll() 280 M*/ 281 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 282 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,0) 283 #else 284 # define TSGLAdaptRegisterDynamic(a,b,c,d) TSGLAdaptRegister(a,b,c,d) 285 #endif 286 287 extern PetscErrorCode TSGLAdaptRegister(const char[],const char[],const char[],PetscErrorCode (*)(TSGLAdapt)); 288 extern PetscErrorCode TSGLAdaptRegisterAll(const char[]); 289 extern PetscErrorCode TSGLAdaptRegisterDestroy(void); 290 extern PetscErrorCode TSGLAdaptInitializePackage(const char[]); 291 extern PetscErrorCode TSGLAdaptFinalizePackage(void); 292 extern PetscErrorCode TSGLAdaptCreate(MPI_Comm,TSGLAdapt*); 293 extern PetscErrorCode TSGLAdaptSetType(TSGLAdapt,const TSGLAdaptType); 294 extern PetscErrorCode TSGLAdaptSetOptionsPrefix(TSGLAdapt,const char[]); 295 extern PetscErrorCode TSGLAdaptChoose(TSGLAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool *); 296 extern PetscErrorCode TSGLAdaptView(TSGLAdapt,PetscViewer); 297 extern PetscErrorCode TSGLAdaptSetFromOptions(TSGLAdapt); 298 extern PetscErrorCode TSGLAdaptDestroy(TSGLAdapt*); 299 300 /*E 301 TSGLAcceptType - String with the name of TSGLAccept scheme or the function 302 with an optional dynamic library name, for example 303 http://www.mcs.anl.gov/petsc/lib.a:mytsglaccept() 304 305 Level: beginner 306 307 .seealso: TSGLSetAcceptType(), TS 308 E*/ 309 #define TSGLAcceptType char* 310 #define TSGLACCEPT_ALWAYS "always" 311 312 typedef PetscErrorCode (*TSGLAcceptFunction)(TS,PetscReal,PetscReal,const PetscReal[],PetscBool *); 313 extern PetscErrorCode TSGLAcceptRegister(const char[],const char[],const char[],TSGLAcceptFunction); 314 315 /*MC 316 TSGLAcceptRegisterDynamic - adds a TSGL acceptance scheme 317 318 Synopsis: 319 PetscErrorCode TSGLAcceptRegisterDynamic(const char *name_scheme,const char *path,const char *name_create,PetscErrorCode (*routine_create)(TS)) 320 321 Not Collective 322 323 Input Parameters: 324 + name_scheme - name of user-defined acceptance scheme 325 . path - path (either absolute or relative) the library containing this scheme 326 . name_create - name of routine to create method context 327 - routine_create - routine to create method context 328 329 Notes: 330 TSGLAcceptRegisterDynamic() may be called multiple times to add several user-defined families. 331 332 If dynamic libraries are used, then the fourth input argument (routine_create) 333 is ignored. 334 335 Sample usage: 336 .vb 337 TSGLAcceptRegisterDynamic("my_scheme",/home/username/my_lib/lib/libO/solaris/mylib.a, 338 "MySchemeCreate",MySchemeCreate); 339 .ve 340 341 Then, your scheme can be chosen with the procedural interface via 342 $ TSGLSetAcceptType(ts,"my_scheme") 343 or at runtime via the option 344 $ -ts_gl_accept_type my_scheme 345 346 Level: advanced 347 348 Notes: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, 349 and others of the form ${any_environmental_variable} occuring in pathname will be 350 replaced with appropriate values. 351 352 .keywords: TSGL, TSGLAcceptType, register 353 354 .seealso: TSGLRegisterAll() 355 M*/ 356 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 357 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,0) 358 #else 359 # define TSGLAcceptRegisterDynamic(a,b,c,d) TSGLAcceptRegister(a,b,c,d) 360 #endif 361 362 /*E 363 TSGLType - family of time integration method within the General Linear class 364 365 Level: beginner 366 367 .seealso: TSGLSetType(), TSGLRegister() 368 E*/ 369 #define TSGLType char* 370 #define TSGL_IRKS "irks" 371 372 /*MC 373 TSGLRegisterDynamic - adds a TSGL implementation 374 375 Synopsis: 376 PetscErrorCode TSGLRegisterDynamic(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 general linear 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 TSGLRegisterDynamic() 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 TSGLRegisterDynamic("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 $ TSGLSetType(ts,"my_scheme") 400 or at runtime via the option 401 $ -ts_gl_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: TSGL, register 410 411 .seealso: TSGLRegisterAll() 412 M*/ 413 #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 414 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,0) 415 #else 416 # define TSGLRegisterDynamic(a,b,c,d) TSGLRegister(a,b,c,d) 417 #endif 418 419 extern PetscErrorCode TSGLRegister(const char[],const char[],const char[],PetscErrorCode(*)(TS)); 420 extern PetscErrorCode TSGLRegisterAll(const char[]); 421 extern PetscErrorCode TSGLRegisterDestroy(void); 422 extern PetscErrorCode TSGLInitializePackage(const char[]); 423 extern PetscErrorCode TSGLFinalizePackage(void); 424 extern PetscErrorCode TSGLSetType(TS,const TSGLType); 425 extern PetscErrorCode TSGLGetAdapt(TS,TSGLAdapt*); 426 extern PetscErrorCode TSGLSetAcceptType(TS,const TSGLAcceptType); 427 428 #define TSARKIMEXType char* 429 #define TSARKIMEX2D "2d" 430 #define TSARKIMEX2E "2e" 431 #define TSARKIMEX3 "3" 432 #define TSARKIMEX4 "4" 433 #define TSARKIMEX5 "5" 434 extern PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType*); 435 extern PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType); 436 extern PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType,PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],const PetscReal[],PetscInt,const PetscReal[],const PetscReal[]); 437 extern PetscErrorCode TSARKIMEXFinalizePackage(void); 438 extern PetscErrorCode TSARKIMEXInitializePackage(const char path[]); 439 extern PetscErrorCode TSARKIMEXRegisterDestroy(void); 440 extern PetscErrorCode TSARKIMEXRegisterAll(void); 441 442 /* 443 PETSc interface to Sundials 444 */ 445 #ifdef PETSC_HAVE_SUNDIALS 446 typedef enum { SUNDIALS_ADAMS=1,SUNDIALS_BDF=2} TSSundialsLmmType; 447 extern const char *TSSundialsLmmTypes[]; 448 typedef enum { SUNDIALS_MODIFIED_GS = 1,SUNDIALS_CLASSICAL_GS = 2 } TSSundialsGramSchmidtType; 449 extern const char *TSSundialsGramSchmidtTypes[]; 450 extern PetscErrorCode TSSundialsSetType(TS,TSSundialsLmmType); 451 extern PetscErrorCode TSSundialsGetPC(TS,PC*); 452 extern PetscErrorCode TSSundialsSetTolerance(TS,PetscReal,PetscReal); 453 extern PetscErrorCode TSSundialsSetMinTimeStep(TS,PetscReal); 454 extern PetscErrorCode TSSundialsSetMaxTimeStep(TS,PetscReal); 455 extern PetscErrorCode TSSundialsGetIterations(TS,PetscInt *,PetscInt *); 456 extern PetscErrorCode TSSundialsSetGramSchmidtType(TS,TSSundialsGramSchmidtType); 457 extern PetscErrorCode TSSundialsSetGMRESRestart(TS,PetscInt); 458 extern PetscErrorCode TSSundialsSetLinearTolerance(TS,PetscReal); 459 extern PetscErrorCode TSSundialsMonitorInternalSteps(TS,PetscBool ); 460 extern PetscErrorCode TSSundialsGetParameters(TS,PetscInt *,long*[],double*[]); 461 #endif 462 463 extern PetscErrorCode TSRKSetTolerance(TS,PetscReal); 464 465 extern PetscErrorCode TSThetaSetTheta(TS,PetscReal); 466 extern PetscErrorCode TSThetaGetTheta(TS,PetscReal*); 467 extern PetscErrorCode TSThetaSetEndpoint(TS,PetscBool); 468 469 extern PetscErrorCode TSAlphaSetAdapt(TS,PetscErrorCode(*)(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*),void*); 470 extern PetscErrorCode TSAlphaAdaptDefault(TS,PetscReal,Vec,Vec,PetscReal*,PetscBool*,void*); 471 extern PetscErrorCode TSAlphaSetRadius(TS,PetscReal); 472 extern PetscErrorCode TSAlphaSetParams(TS,PetscReal,PetscReal,PetscReal); 473 extern PetscErrorCode TSAlphaGetParams(TS,PetscReal*,PetscReal*,PetscReal*); 474 475 extern PetscErrorCode TSSetDM(TS,DM); 476 extern PetscErrorCode TSGetDM(TS,DM*); 477 478 extern PetscErrorCode SNESTSFormFunction(SNES,Vec,Vec,void*); 479 extern PetscErrorCode SNESTSFormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*); 480 481 PETSC_EXTERN_CXX_END 482 #endif 483