Actual source code: tsimpl.h

  1:   /* $Id: tsimpl.h,v 1.25 2001/09/07 20:12:01 bsmith Exp $ */

  3: #ifndef __TSIMPL_H

 6:  #include petscts.h

  8: /*
  9:     Timesteping context. 
 10:       General case: U_t = F(t,U) <-- the right-hand-side function
 11:       Linear  case: U_t = A(t) U <-- the right-hand-side matrix
 12:       Linear (no time) case: U_t = A U <-- the right-hand-side matrix
 13: */

 15: /*
 16:      Maximum number of monitors you can run with a single TS
 17: */
 18: #define MAXTSMONITORS 5 

 20: struct _TSOps {
 21:   int (*rhsmatrix)(TS, PetscReal, Mat *, Mat *, MatStructure *, void *),
 22:       (*rhsfunction)(TS, PetscReal, Vec, Vec, void *),
 23:       (*rhsjacobian)(TS, PetscReal, Vec, Mat *, Mat *, MatStructure *, void *),
 24:       (*applymatrixbc)(TS, Mat, Mat, void *),
 25:       (*rhsbc)(TS, PetscReal, Vec, void *),
 26:       (*applyrhsbc)(TS, Vec, void *),
 27:       (*applysolbc)(TS, Vec, void *),
 28:       (*prestep)(TS),
 29:       (*update)(TS, PetscReal, PetscReal *),
 30:       (*poststep)(TS),
 31:       (*reform)(TS),
 32:       (*reallocate)(TS),
 33:       (*setup)(TS),
 34:       (*step)(TS,int *, PetscReal *),
 35:       (*setfromoptions)(TS),
 36:       (*printhelp)(TS, char *),
 37:       (*destroy)(TS),
 38:       (*view)(TS, PetscViewer);
 39: };

 41: struct _p_TS {
 42:   PETSCHEADER(struct _TSOps)
 43:   TSProblemType problem_type;
 44:   Vec           vec_sol, vec_sol_always;

 46:   /* ---------------- User (or PETSc) Provided stuff ---------------------*/
 47:   int  (*monitor[MAXTSMONITORS])(TS,int,PetscReal,Vec,void*); /* returns control to user after */
 48:   int  (*mdestroy[MAXTSMONITORS])(void*);
 49:   void *monitorcontext[MAXTSMONITORS];                 /* residual calculation, allows user */
 50:   int  numbermonitors;                                 /* to, for instance, print residual norm, etc. */

 52:   /* Identifies this as a grid TS structure */
 53:   PetscTruth  isGTS;                                 /* This problem arises from an underlying grid */
 54:   PetscTruth *isExplicit;                            /* Indicates which fields have explicit time dependence */
 55:   int        *Iindex;                                /* The index of the identity for each time dependent field */

 57:   /* ---------------------Linear Iteration---------------------------------*/
 58:   SLES sles;
 59:   Mat  A, B;                                         /* user provided matrix and preconditioner */

 61:   /* ---------------------Nonlinear Iteration------------------------------*/
 62:   SNES  snes;
 63:   void *funP;
 64:   void *jacP;
 65:   void *bcP;


 68:   /* --- Data that is unique to each particular solver --- */
 69:   int   setupcalled;            /* true if setup has been called */
 70:   void *data;                   /* implementationspecific data */
 71:   void *user;                   /* user context */

 73:   /* ------------------  Parameters -------------------------------------- */
 74:   int       max_steps;              /* max number of steps */
 75:   PetscReal max_time;               /* max time allowed */
 76:   PetscReal time_step;              /* current time increment */
 77:   PetscReal time_step_old;          /* previous time increment */
 78:   PetscReal initial_time_step;      /* initial time increment */
 79:   int       steps;                  /* steps taken so far */
 80:   PetscReal ptime;                  /* time taken so far */
 81:   int       linear_its;             /* total number of linear solver iterations */
 82:   int       nonlinear_its;          /* total number of nonlinear solver iterations */

 84:   /* ------------------- Default work-area management ------------------ */
 85:   int  nwork;
 86:   Vec *work;
 87: };

 89: EXTERN int TSMonitor(TS,int,PetscReal,Vec);
 90: EXTERN int TSComputeRHSBoundaryConditions(TS,PetscReal,Vec);

 92: #endif