Actual source code: snesimpl.h
petsc-3.9.4 2018-09-11
2: #ifndef __SNESIMPL_H
5: #include <petscsnes.h>
6: #include <petsc/private/petscimpl.h>
8: PETSC_EXTERN PetscBool SNESRegisterAllCalled;
9: PETSC_EXTERN PetscErrorCode SNESRegisterAll(void);
11: typedef struct _SNESOps *SNESOps;
13: struct _SNESOps {
14: PetscErrorCode (*computeinitialguess)(SNES,Vec,void*);
15: PetscErrorCode (*computescaling)(Vec,Vec,void*);
16: PetscErrorCode (*update)(SNES, PetscInt); /* General purpose function for update */
17: PetscErrorCode (*converged)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
18: PetscErrorCode (*convergeddestroy)(void*);
19: PetscErrorCode (*setup)(SNES); /* routine to set up the nonlinear solver */
20: PetscErrorCode (*solve)(SNES); /* actual nonlinear solver */
21: PetscErrorCode (*view)(SNES,PetscViewer);
22: PetscErrorCode (*setfromoptions)(PetscOptionItems*,SNES); /* sets options from database */
23: PetscErrorCode (*destroy)(SNES);
24: PetscErrorCode (*reset)(SNES);
25: PetscErrorCode (*usercompute)(SNES,void**);
26: PetscErrorCode (*userdestroy)(void**);
27: PetscErrorCode (*computevariablebounds)(SNES,Vec,Vec); /* user provided routine to set box constrained variable bounds */
28: PetscErrorCode (*computepfunction)(SNES,Vec,Vec,void*);
29: PetscErrorCode (*computepjacobian)(SNES,Vec,Mat,Mat,void*);
30: PetscErrorCode (*load)(SNES,PetscViewer);
31: };
33: /*
34: Nonlinear solver context
35: */
36: #define MAXSNESMONITORS 5
38: struct _p_SNES {
39: PETSCHEADER(struct _SNESOps);
40: DM dm;
41: PetscBool dmAuto; /* SNES created currently used DM automatically */
42: SNES npc;
43: PCSide npcside;
44: PetscBool usesnpc; /* type can use a nonlinear preconditioner */
46: /* ------------------------ User-provided stuff -------------------------------*/
47: void *user; /* user-defined context */
49: Vec vec_rhs; /* If non-null, solve F(x) = rhs */
50: Vec vec_sol; /* pointer to solution */
52: Vec vec_func; /* pointer to function */
54: Mat jacobian; /* Jacobian matrix */
55: Mat jacobian_pre; /* preconditioner matrix */
56: void *initialguessP; /* user-defined initial guess context */
57: KSP ksp; /* linear solver context */
58: SNESLineSearch linesearch; /* line search context */
59: PetscBool usesksp;
60: MatStructure matstruct; /* Used by Picard solver */
62: Vec vec_sol_update; /* pointer to solution update */
64: Vec scaling; /* scaling vector */
65: void *scaP; /* scaling context */
67: PetscReal precheck_picard_angle; /* For use with SNESLineSearchPreCheckPicard */
69: /* ------------------------Time stepping hooks-----------------------------------*/
71: /* ---------------- PETSc-provided (or user-provided) stuff ---------------------*/
73: PetscErrorCode (*monitor[MAXSNESMONITORS])(SNES,PetscInt,PetscReal,void*); /* monitor routine */
74: PetscErrorCode (*monitordestroy[MAXSNESMONITORS])(void**); /* monitor context destroy routine */
75: void *monitorcontext[MAXSNESMONITORS]; /* monitor context */
76: PetscInt numbermonitors; /* number of monitors */
77: void *cnvP; /* convergence context */
78: SNESConvergedReason reason;
79: PetscBool errorifnotconverged;
81: /* --- Routines and data that are unique to each particular solver --- */
83: PetscBool setupcalled; /* true if setup has been called */
84: void *data; /* implementation-specific data */
86: /* -------------------------- Parameters -------------------------------------- */
88: PetscInt max_its; /* max number of iterations */
89: PetscInt max_funcs; /* max number of function evals */
90: PetscInt nfuncs; /* number of function evaluations */
91: PetscInt iter; /* global iteration number */
92: PetscInt linear_its; /* total number of linear solver iterations */
93: PetscReal norm; /* residual norm of current iterate */
94: PetscReal rtol; /* relative tolerance */
95: PetscReal divtol; /* relative divergence tolerance */
96: PetscReal abstol; /* absolute tolerance */
97: PetscReal stol; /* step length tolerance*/
98: PetscReal deltatol; /* trust region convergence tolerance */
99: PetscBool forceiteration; /* Force SNES to take at least one iteration regardless of the initial residual norm */
100: PetscInt lagpreconditioner; /* SNESSetLagPreconditioner() */
101: PetscInt lagjacobian; /* SNESSetLagJacobian() */
102: PetscInt jac_iter; /* The present iteration of the Jacobian lagging */
103: PetscBool lagjac_persist; /* The jac_iter persists until reset */
104: PetscInt pre_iter; /* The present iteration of the Preconditioner lagging */
105: PetscBool lagpre_persist; /* The pre_iter persists until reset */
106: PetscInt gridsequence; /* number of grid sequence steps to take; defaults to zero */
108: PetscBool tolerancesset; /* SNESSetTolerances() called and tolerances should persist through SNESCreate_XXX()*/
110: PetscBool vec_func_init_set; /* the initial function has been set */
112: SNESNormSchedule normschedule; /* Norm computation type for SNES instance */
113: SNESFunctionType functype; /* Function type for the SNES instance */
115: /* ------------------------ Default work-area management ---------------------- */
117: PetscInt nwork;
118: Vec *work;
120: /* ------------------------- Miscellaneous Information ------------------------ */
122: PetscReal *conv_hist; /* If !0, stores function norm (or
123: gradient norm) at each iteration */
124: PetscInt *conv_hist_its; /* linear iterations for each Newton step */
125: PetscInt conv_hist_len; /* size of convergence history array */
126: PetscInt conv_hist_max; /* actual amount of data in conv_history */
127: PetscBool conv_hist_reset; /* reset counter for each new SNES solve */
128: PetscBool conv_malloc;
130: PetscBool counters_reset; /* reset counter for each new SNES solve */
132: /* the next two are used for failures in the line search; they should be put elsewhere */
133: PetscInt numFailures; /* number of unsuccessful step attempts */
134: PetscInt maxFailures; /* maximum number of unsuccessful step attempts */
136: PetscInt numLinearSolveFailures;
137: PetscInt maxLinearSolveFailures;
139: PetscBool domainerror; /* set with SNESSetFunctionDomainError() */
141: PetscBool ksp_ewconv; /* flag indicating use of Eisenstat-Walker KSP convergence criteria */
142: void *kspconvctx; /* Eisenstat-Walker KSP convergence context */
144: /* SNESConvergedDefault context: split it off into a separate var/struct to be passed as context to SNESConvergedDefault? */
145: PetscReal ttol; /* rtol*initial_residual_norm */
146: PetscReal rnorm0; /* initial residual norm (used for divergence testing) */
147:
148: Vec *vwork; /* more work vectors for Jacobian approx */
149: PetscInt nvwork;
151: PetscBool mf; /* -snes_mf was used on this snes */
152: PetscBool mf_operator; /* -snes_mf_operator was used on this snes */
153: PetscInt mf_version; /* The version of snes_mf used */
155: PetscReal vizerotolerance; /* tolerance for considering an x[] value to be on the bound */
156: Vec xl,xu; /* upper and lower bounds for box constrained VI problems */
157: PetscInt ntruebounds; /* number of non-infinite bounds set for VI box constraints */
158: PetscBool usersetbounds; /* bounds have been set via SNESVISetVariableBounds(), rather than via computevariablebounds() callback. */
160: PetscBool alwayscomputesfinalresidual; /* Does SNESSolve_XXX always compute the value of the residual at the final
161: * solution and put it in vec_func? Used inside SNESSolve_FAS to determine
162: * if the final residual must be computed before restricting or prolonging
163: * it. */
165: };
167: typedef struct _p_DMSNES *DMSNES;
168: typedef struct _DMSNESOps *DMSNESOps;
169: struct _DMSNESOps {
170: PetscErrorCode (*computefunction)(SNES,Vec,Vec,void*);
171: PetscErrorCode (*computejacobian)(SNES,Vec,Mat,Mat,void*);
173: /* objective */
174: PetscErrorCode (*computeobjective)(SNES,Vec,PetscReal*,void*);
176: /* Picard iteration functions */
177: PetscErrorCode (*computepfunction)(SNES,Vec,Vec,void*);
178: PetscErrorCode (*computepjacobian)(SNES,Vec,Mat,Mat,void*);
180: /* User-defined smoother */
181: PetscErrorCode (*computegs)(SNES,Vec,Vec,void*);
183: PetscErrorCode (*destroy)(DMSNES);
184: PetscErrorCode (*duplicate)(DMSNES,DMSNES);
185: };
187: struct _p_DMSNES {
188: PETSCHEADER(struct _DMSNESOps);
189: void *functionctx;
190: void *gsctx;
191: void *pctx;
192: void *jacobianctx;
193: void *objectivectx;
195: void *data;
197: /* This is NOT reference counted. The DM on which this context was first created is cached here to implement one-way
198: * copy-on-write. When DMGetDMSNESWrite() sees a request using a different DM, it makes a copy. Thus, if a user
199: * only interacts directly with one level, e.g., using SNESSetFunction(), then SNESSetUp_FAS() is called to build
200: * coarse levels, then the user changes the routine with another call to SNESSetFunction(), it automatically
201: * propagates to all the levels. If instead, they get out a specific level and set the function on that level,
202: * subsequent changes to the original level will no longer propagate to that level.
203: */
204: DM originaldm;
205: };
206: PETSC_EXTERN PetscErrorCode DMGetDMSNES(DM,DMSNES*);
207: PETSC_EXTERN PetscErrorCode DMSNESView(DMSNES,PetscViewer);
208: PETSC_EXTERN PetscErrorCode DMSNESLoad(DMSNES,PetscViewer);
209: PETSC_EXTERN PetscErrorCode DMGetDMSNESWrite(DM,DMSNES*);
212: /* Context for Eisenstat-Walker convergence criteria for KSP solvers */
213: typedef struct {
214: PetscInt version; /* flag indicating version 1 or 2 of test */
215: PetscReal rtol_0; /* initial rtol */
216: PetscReal rtol_last; /* last rtol */
217: PetscReal rtol_max; /* maximum rtol */
218: PetscReal gamma; /* mult. factor for version 2 rtol computation */
219: PetscReal alpha; /* power for version 2 rtol computation */
220: PetscReal alpha2; /* power for safeguard */
221: PetscReal threshold; /* threshold for imposing safeguard */
222: PetscReal lresid_last; /* linear residual from last iteration */
223: PetscReal norm_last; /* function norm from last iteration */
224: PetscReal norm_first; /* function norm from the beginning of the first iteration. */
225: } SNESKSPEW;
227: PETSC_STATIC_INLINE PetscErrorCode SNESLogConvergenceHistory(SNES snes,PetscReal res,PetscInt its)
228: {
232: PetscObjectSAWsTakeAccess((PetscObject)snes);
233: if (snes->conv_hist && snes->conv_hist_max > snes->conv_hist_len) {
234: if (snes->conv_hist) snes->conv_hist[snes->conv_hist_len] = res;
235: if (snes->conv_hist_its) snes->conv_hist_its[snes->conv_hist_len] = its;
236: snes->conv_hist_len++;
237: }
238: PetscObjectSAWsGrantAccess((PetscObject)snes);
239: return(0);
240: }
242: PETSC_EXTERN PetscErrorCode SNESVIProjectOntoBounds(SNES,Vec);
243: PETSC_INTERN PetscErrorCode SNESVICheckLocalMin_Private(SNES,Mat,Vec,Vec,PetscReal,PetscBool*);
244: PETSC_INTERN PetscErrorCode SNESReset_VI(SNES);
245: PETSC_INTERN PetscErrorCode SNESDestroy_VI(SNES);
246: PETSC_INTERN PetscErrorCode SNESView_VI(SNES,PetscViewer);
247: PETSC_INTERN PetscErrorCode SNESSetFromOptions_VI(PetscOptionItems*,SNES);
248: PETSC_INTERN PetscErrorCode SNESSetUp_VI(SNES);
249: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*SNESVIComputeVariableBoundsFunction)(SNES,Vec,Vec);
250: PETSC_INTERN PetscErrorCode SNESVISetComputeVariableBounds_VI(SNES,SNESVIComputeVariableBoundsFunction);
251: PETSC_INTERN PetscErrorCode SNESVISetVariableBounds_VI(SNES,Vec,Vec);
252: PETSC_INTERN PetscErrorCode SNESConvergedDefault_VI(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*);
254: PetscErrorCode SNESScaleStep_Private(SNES,Vec,PetscReal*,PetscReal*,PetscReal*,PetscReal*);
255: PETSC_EXTERN PetscErrorCode DMSNESCheckFromOptions_Internal(SNES,DM,Vec,Vec,PetscErrorCode (**)(PetscInt,PetscReal,const PetscReal[],PetscInt,PetscScalar*,void*),void**);
257: PETSC_EXTERN PetscLogEvent SNES_Solve;
258: PETSC_EXTERN PetscLogEvent SNES_LineSearch;
259: PETSC_EXTERN PetscLogEvent SNES_FunctionEval;
260: PETSC_EXTERN PetscLogEvent SNES_JacobianEval;
261: PETSC_EXTERN PetscLogEvent SNES_NGSEval;
262: PETSC_EXTERN PetscLogEvent SNES_NGSFuncEval;
263: PETSC_EXTERN PetscLogEvent SNES_NPCSolve;
264: PETSC_EXTERN PetscLogEvent SNES_ObjectiveEval;
266: extern PetscBool SNEScite;
267: extern const char SNESCitation[];
269: /*
270: Either generate an error or mark as diverged when a real from a SNES function norm is Nan or Inf
271: */
272: #define SNESCheckFunctionNorm(snes,beta) \
273: if (PetscIsInfOrNanReal(beta)) {\
274: if (snes->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_NOT_CONVERGED,"SNESSolve has not converged due to Nan or Inf norm");\
275: else {\
276: PetscBool domainerror;\
277: PetscErrorCode MPIU_Allreduce((int*)&snes->domainerror,(int*)&domainerror,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)snes));\
278: if (domainerror) snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;\
279: else snes->reason = SNES_DIVERGED_FNORM_NAN;\
280: return(0);\
281: }\
282: }
284: #define SNESCheckKSPSolve(snes)\
285: {\
286: KSPConvergedReason kspreason; \
287: PetscErrorCode KSPGetConvergedReason(snes->ksp,&kspreason);\
288: if (kspreason < 0) {\
289: if (kspreason == KSP_DIVERGED_NANORINF) {\
290: PetscBool domainerror;\
291: MPIU_Allreduce((int*)&snes->domainerror,(int*)&domainerror,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)snes)); \
292: if (domainerror) snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;\
293: else snes->reason = SNES_DIVERGED_LINEAR_SOLVE; \
294: return(0);\
295: } else {\
296: if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {\
297: PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);\
298: snes->reason = SNES_DIVERGED_LINEAR_SOLVE;\
299: return(0);\
300: }\
301: }\
302: }\
303: }
305: #endif