Actual source code: snesimpl.h
2: #ifndef __SNESIMPL_H
5: #include petscsnes.h
7: typedef struct _SNESOps *SNESOps;
9: struct _SNESOps {
10: PetscErrorCode (*computefunction)(SNES,Vec,Vec,void*);
11: PetscErrorCode (*computejacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
12: PetscErrorCode (*computescaling)(Vec,Vec,void*);
13: PetscErrorCode (*update)(SNES, PetscInt); /* General purpose function for update */
14: PetscErrorCode (*converged)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
15: PetscErrorCode (*setup)(SNES); /* routine to set up the nonlinear solver */
16: PetscErrorCode (*solve)(SNES); /* actual nonlinear solver */
17: PetscErrorCode (*view)(SNES,PetscViewer);
18: PetscErrorCode (*setfromoptions)(SNES); /* sets options from database */
19: PetscErrorCode (*destroy)(SNES);
20: };
22: /*
23: Nonlinear solver context
24: */
25: #define MAXSNESMONITORS 5
27: struct _p_SNES {
28: PETSCHEADER(struct _SNESOps);
30: /* ------------------------ User-provided stuff -------------------------------*/
31: void *user; /* user-defined context */
33: Vec vec_sol,vec_sol_always; /* pointer to solution */
34: Vec vec_sol_update_always; /* pointer to solution update */
36: Vec vec_func,vec_func_always; /* pointer to function */
37: Vec afine; /* If non-null solve F(x) = afine */
38: void *funP; /* user-defined function context */
41: Mat jacobian; /* Jacobian matrix */
42: Mat jacobian_pre; /* preconditioner matrix */
43: void *jacP; /* user-defined Jacobian context */
44: KSP ksp; /* linear solver context */
46: Vec scaling; /* scaling vector */
47: void *scaP; /* scaling context */
49: /* ------------------------Time stepping hooks-----------------------------------*/
51: /* ---------------- PETSc-provided (or user-provided) stuff ---------------------*/
53: PetscErrorCode (*monitor[MAXSNESMONITORS])(SNES,PetscInt,PetscReal,void*); /* monitor routine */
54: PetscErrorCode (*monitordestroy[MAXSNESMONITORS])(void*); /* monitor context destroy routine */
55: void *monitorcontext[MAXSNESMONITORS]; /* monitor context */
56: PetscInt numbermonitors; /* number of monitors */
57: void *cnvP; /* convergence context */
58: SNESConvergedReason reason;
60: /* --- Routines and data that are unique to each particular solver --- */
62: PetscTruth setupcalled; /* true if setup has been called */
63: void *data; /* implementation-specific data */
65: /* -------------------------- Parameters -------------------------------------- */
67: PetscInt max_its; /* max number of iterations */
68: PetscInt max_funcs; /* max number of function evals */
69: PetscInt nfuncs; /* number of function evaluations */
70: PetscInt iter; /* global iteration number */
71: PetscInt linear_its; /* total number of linear solver iterations */
72: PetscReal norm; /* residual norm of current iterate */
73: PetscReal rtol; /* relative tolerance */
74: PetscReal abstol; /* absolute tolerance */
75: PetscReal xtol; /* relative tolerance in solution */
76: PetscReal deltatol; /* trust region convergence tolerance */
77: PetscTruth printreason; /* print reason for convergence/divergence after each solve */
78: /* ------------------------ Default work-area management ---------------------- */
80: PetscInt nwork;
81: Vec *work;
83: /* ------------------------- Miscellaneous Information ------------------------ */
85: PetscReal *conv_hist; /* If !0, stores function norm (or
86: gradient norm) at each iteration */
87: PetscInt *conv_hist_its; /* linear iterations for each Newton step */
88: PetscInt conv_hist_len; /* size of convergence history array */
89: PetscInt conv_hist_max; /* actual amount of data in conv_history */
90: PetscTruth conv_hist_reset; /* reset counter for each new SNES solve */
92: /* the next two are used for failures in the line search; they should be put into the LS struct */
93: PetscInt numFailures; /* number of unsuccessful step attempts */
94: PetscInt maxFailures; /* maximum number of unsuccessful step attempts */
96: PetscInt numLinearSolveFailures;
97: PetscInt maxLinearSolveFailures;
99: /*
100: These are REALLY ugly and don't belong here, but since they must
101: be destroyed at the conclusion we have to put them somewhere.
102: */
103: PetscTruth ksp_ewconv; /* flag indicating use of Eisenstat-Walker KSP convergence criteria */
104: void *kspconvctx; /* KSP convergence context */
106: PetscReal ttol; /* used by default convergence test routine */
108: Vec *vwork; /* more work vectors for Jacobian approx */
109: PetscInt nvwork;
110: };
112: /* Context for Eisenstat-Walker convergence criteria for KSP solvers */
113: typedef struct {
114: PetscInt version; /* flag indicating version 1 or 2 of test */
115: PetscReal rtol_0; /* initial rtol */
116: PetscReal rtol_last; /* last rtol */
117: PetscReal rtol_max; /* maximum rtol */
118: PetscReal gamma; /* mult. factor for version 2 rtol computation */
119: PetscReal alpha; /* power for version 2 rtol computation */
120: PetscReal alpha2; /* power for safeguard */
121: PetscReal threshold; /* threshold for imposing safeguard */
122: PetscReal lresid_last; /* linear residual from last iteration */
123: PetscReal norm_last; /* function norm from last iteration */
124: } SNESKSPEW;
126: #define SNESLogConvHistory(snes,res,its) \
127: { if (snes->conv_hist && snes->conv_hist_max > snes->conv_hist_len) \
128: { if (snes->conv_hist) snes->conv_hist[snes->conv_hist_len] = res; \
129: if (snes->conv_hist_its) snes->conv_hist_its[snes->conv_hist_len] = its; \
130: snes->conv_hist_len++;\
131: }}
133: #define SNESMonitor(snes,it,rnorm) \
134: { PetscErrorCode _ierr; PetscInt _i,_im = snes->numbermonitors; \
135: for (_i=0; _i<_im; _i++) {\
136: _(*snes->monitor[_i])(snes,it,rnorm,snes->monitorcontext[_i]);CHKERRQ(_ierr); \
137: } \
138: }
140: PetscErrorCode SNES_KSPSolve(SNES,KSP,Vec,Vec);
141: PetscErrorCode SNESScaleStep_Private(SNES,Vec,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
145: #endif