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