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