Actual source code: tsimpl.h

  2: #ifndef __TSIMPL_H

  5: #include <petscts.h>

  7: /*
  8:     Timesteping context.
  9:       General DAE: F(t,U,U_t) = 0, required Jacobian is G'(U) where G(U) = F(t,U,U0+a*U)
 10:       General ODE: U_t = F(t,U) <-- the right-hand-side function
 11:       Linear  ODE: U_t = A(t) U <-- the right-hand-side matrix
 12:       Linear (no time) ODE: 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: typedef struct _TSOps *TSOps;

 22: struct _TSOps {
 23:   PetscErrorCode (*snesfunction)(SNES,Vec,Vec,TS);
 24:   PetscErrorCode (*snesjacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,TS);
 25:   PetscErrorCode (*prestep)(TS);
 26:   PetscErrorCode (*poststep)(TS);
 27:   PetscErrorCode (*setup)(TS);
 28:   PetscErrorCode (*step)(TS);
 29:   PetscErrorCode (*solve)(TS);
 30:   PetscErrorCode (*interpolate)(TS,PetscReal,Vec);
 31:   PetscErrorCode (*setfromoptions)(TS);
 32:   PetscErrorCode (*destroy)(TS);
 33:   PetscErrorCode (*view)(TS,PetscViewer);
 34:   PetscErrorCode (*reset)(TS);
 35: };

 37: struct _TSUserOps {
 38:   PetscErrorCode (*rhsfunction)(TS,PetscReal,Vec,Vec,void*);
 39:   PetscErrorCode (*rhsjacobian)(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
 40:   PetscErrorCode (*ifunction)(TS,PetscReal,Vec,Vec,Vec,void*);
 41:   PetscErrorCode (*ijacobian)(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
 42: };

 44: struct _p_TS {
 45:   PETSCHEADER(struct _TSOps);

 47:   struct _TSUserOps *userops;
 48:   DM            dm;
 49:   TSProblemType problem_type;
 50:   Vec           vec_sol;

 52:   /* ---------------- User (or PETSc) Provided stuff ---------------------*/
 53:   PetscErrorCode (*monitor[MAXTSMONITORS])(TS,PetscInt,PetscReal,Vec,void*); /* returns control to user after */
 54:   PetscErrorCode (*mdestroy[MAXTSMONITORS])(void**);
 55:   void *monitorcontext[MAXTSMONITORS];                 /* residual calculation, allows user */
 56:   PetscInt  numbermonitors;                                 /* to, for instance, print residual norm, etc. */

 58:   /* ---------------------- IMEX support ---------------------------------*/
 59:   /* These extra slots are only used when the user provides both Implicit and RHS */
 60:   Mat Arhs;     /* Right hand side matrix */
 61:   Mat Brhs;     /* Right hand side preconditioning matrix */
 62:   Vec Frhs;     /* Right hand side function value */

 64:   /* This is a general caching scheme to avoid recomputing the Jacobian at a place that has been previously been evaluated.
 65:    * The present use case is that TSComputeRHSFunctionLinear() evaluates the Jacobian once and we don't want it to be immeditely re-evaluated.
 66:    */
 67:   struct {
 68:     PetscReal time;             /* The time at which the matrices were last evaluated */
 69:     Vec X;                      /* Solution vector at which the Jacobian was last evaluated */
 70:     PetscInt Xstate;            /* State of the solution vector */
 71:     MatStructure mstructure;    /* The structure returned */
 72:   } rhsjacobian;

 74:   struct {
 75:     PetscReal time;             /* The time at which the matrices were last evaluated */
 76:     Vec X;                      /* Solution vector at which the Jacobian was last evaluated */
 77:     Vec Xdot;                   /* Time derivative of the state vector at which the Jacobian was last evaluated */
 78:     PetscInt Xstate;            /* State of the solution vector */
 79:     PetscInt Xdotstate;         /* State of the solution vector */
 80:     MatStructure mstructure;    /* The structure returned */
 81:     PetscReal shift;            /* The derivative of the lhs wrt to Xdot */
 82:     PetscBool imex;             /* Flag of the method if it was started as an imex method */
 83:   } ijacobian;

 85:   /* ---------------------Nonlinear Iteration------------------------------*/
 86:   SNES  snes;
 87:   void *funP;
 88:   void *jacP;


 91:   /* --- Data that is unique to each particular solver --- */
 92:   PetscInt setupcalled;             /* true if setup has been called */
 93:   void     *data;                   /* implementationspecific data */
 94:   void     *user;                   /* user context */

 96:   /* ------------------  Parameters -------------------------------------- */
 97:   PetscInt  max_steps;              /* max number of steps */
 98:   PetscReal max_time;               /* max time allowed */
 99:   PetscReal time_step;              /* current/completed time increment */
100:   PetscReal time_step_prev;         /* previous time step  */
101:   PetscInt  steps;                  /* steps taken so far */
102:   PetscReal ptime;                  /* time at the start of the current step (stage time is internal if it exists) */
103:   PetscInt  linear_its;             /* total number of linear solver iterations */
104:   PetscInt  nonlinear_its;          /* total number of nonlinear solver iterations */

106:   PetscInt num_snes_failures;
107:   PetscInt max_snes_failures;
108:   TSConvergedReason reason;
109:   PetscBool errorifstepfailed;
110:   PetscInt  exact_final_time;   /* PETSC_DECIDE, PETSC_TRUE, or PETSC_FALSE */
111:   PetscBool retain_stages;
112:   PetscInt reject,max_reject;

114:   /* ------------------- Default work-area management ------------------ */
115:   PetscInt nwork;
116:   Vec      *work;
117: };



122: #endif