xref: /petsc/include/petsc/private/linesearchimpl.h (revision bd89dbf26d8a5efecb980364933175da61864cd7) !
1 #pragma once
2 
3 #include <petscsnes.h>
4 #include <petsc/private/petscimpl.h>
5 
6 PETSC_EXTERN PetscBool      SNESLineSearchRegisterAllCalled;
7 PETSC_EXTERN PetscErrorCode SNESLineSearchRegisterAll(void);
8 PETSC_EXTERN PetscLogEvent  SNESLINESEARCH_Apply;
9 
10 typedef struct _LineSearchOps *LineSearchOps;
11 
12 struct _LineSearchOps {
13   PetscErrorCode (*view)(SNESLineSearch, PetscViewer);
14   SNESLineSearchApplyFn *apply;
15   PetscErrorCode (*precheck)(SNESLineSearch, Vec, Vec, PetscBool *, void *);
16   SNESLineSearchVIProjectFn  *viproject;
17   SNESLineSearchVINormFn     *vinorm;
18   SNESLineSearchVIDirDerivFn *vidirderiv;
19   PetscErrorCode (*postcheck)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
20   PetscErrorCode (*setfromoptions)(SNESLineSearch, PetscOptionItems);
21   PetscErrorCode (*reset)(SNESLineSearch);
22   PetscErrorCode (*destroy)(SNESLineSearch);
23   PetscErrorCode (*setup)(SNESLineSearch);
24   PetscErrorCode (*snesfunc)(SNES, Vec, Vec);
25 };
26 
27 #define MAXSNESLSMONITORS 5
28 
29 struct _p_LineSearch {
30   PETSCHEADER(struct _LineSearchOps);
31 
32   SNES snes;
33 
34   void *data;
35 
36   PetscBool setupcalled;
37 
38   Vec vec_sol;
39   Vec vec_sol_new;
40   Vec vec_func;
41   Vec vec_func_new;
42   Vec vec_update;
43 
44   PetscInt nwork;
45   Vec     *work;
46 
47   PetscReal lambda;
48 
49   PetscBool norms;
50   PetscReal fnorm;
51   PetscReal ynorm;
52   PetscReal xnorm;
53   PetscBool keeplambda;
54 
55   PetscReal damping;
56   PetscReal maxlambda;
57   PetscReal minlambda;
58   PetscInt  max_it;
59   PetscReal rtol;
60   PetscReal atol;
61   PetscReal ltol;
62   PetscInt  order;
63 
64   PetscReal precheck_picard_angle;
65 
66   void *precheckctx;
67   void *postcheckctx;
68 
69   PetscBool jacobiandomainerror; /* set with SNESSetJacobianDomainError() */
70   PetscBool checkjacdomainerror; /* does it check Jacobian domain error after Jacobian evaluations */
71 
72   SNESLineSearchReason reason;
73 
74   PetscViewer monitor;
75   PetscErrorCode (*monitorftns[MAXSNESLSMONITORS])(SNESLineSearch, void *); /* monitor routine */
76   PetscCtxDestroyFn *monitordestroy[MAXSNESLSMONITORS];                     /* monitor context destroy routine */
77   void              *monitorcontext[MAXSNESLSMONITORS];                     /* monitor context */
78   PetscInt           numbermonitors;                                        /* number of monitors */
79 };
80 
81 /*MC
82   SNESLineSearchCheckFunctionDomainError - Called after a `SNESComputeFunction()` and `VecNorm()` in a `SNES` line search to check if the function norm is infinity or NaN and
83   if the function callback set with `SNESSetFunction()` called `SNESSetFunctionDomainError()`.
84 
85   Synopsis:
86   #include <snesimpl.h>
87   void SNESLineSearchCheckFunctionDomainError(SNES snes, SNESLineSearch ls, PetscReal fnorm)
88 
89   Collective
90 
91   Input Parameters:
92 +  snes  - the `SNES` object
93 .  ls    - the `SNESLineSearch` object
94 -  fnorm - the value of the norm
95 
96  Level: developer
97 
98  Notes:
99  If `fnorm` is infinity or NaN and `SNESSetErrorIfNotConverged()` was set, this immediately generates a `PETSC_ERR_CONV_FAILED`.
100 
101  If `fnorm` is infinity or NaN and `SNESSetFunctionDomainError()` was called, this sets the `SNESLineSearchReason` to `SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN`
102  and exits the solver
103 
104  Otherwise, if `fnorm` is infinity or NaN, this sets the `SNESLineSearchReason` to `SNES_LINESEARCH_FAILED_NANORINF` and exits the line search
105 
106  See `SNESCheckFunctionDomainError()` for an explanation of the design
107 
108 .seealso: [](ch_snes), `SNESCheckFunctionDomainError()`, `SNESSetFunctionDomainError()`, `PETSC_ERR_CONV_FAILED`, `SNESSetErrorIfNotConverged()`, `SNES_DIVERGED_FUNCTION_DOMAIN`,
109           `SNESConvergedReason`, `SNES_DIVERGED_FUNCTION_NAN`
110 MC*/
111 #define SNESLineSearchCheckFunctionDomainError(snes, ls, fnorm) \
112   do { \
113     if (PetscIsInfOrNanReal(fnorm)) { \
114       PetscCheck(!snes->errorifnotconverged, PetscObjectComm((PetscObject)ls), PETSC_ERR_NOT_CONVERGED, "SNES line search failure due to infinity or NaN norm"); \
115       { \
116         PetscBool functiondomainerror; \
117         PetscCallMPI(MPIU_Allreduce(&snes->functiondomainerror, &functiondomainerror, 1, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)ls))); \
118         if (functiondomainerror) { \
119           ls->reason                = SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN; \
120           snes->functiondomainerror = PETSC_FALSE; \
121         } else ls->reason = SNES_LINESEARCH_FAILED_NANORINF; \
122         PetscFunctionReturn(PETSC_SUCCESS); \
123       } \
124     } \
125   } while (0)
126 
127 /*MC
128   SNESLineSearchCheckObjectiveDomainError - Called after a `SNESComputeObjective()` in a `SNESLineSearch` to check if the objective value is infinity or NaN and/or
129   if the function callback set with `SNESSetObjective()` called `SNESSetObjectiveDomainError()`.
130 
131   Synopsis:
132   #include <snesimpl.h>
133   void SNESLineSearchCheckObjectiveDomainError(SNES snes, PetscReal fobj)
134 
135   Collective
136 
137   Input Parameters:
138 + snes - the `SNES` solver object
139 - fobj - the value of the objective function
140 
141   Level: developer
142 
143   Notes:
144   If `fobj` is infinity or NaN and `SNESSetErrorIfNotConverged()` was set, this immediately generates a `PETSC_ERR_CONV_FAILED`.
145 
146   If `SNESSetObjectiveDomainError()` was called, this sets the `SNESLineSearchReason` to `SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN`
147   and exits the line search
148 
149   Otherwise if `fobj` is infinity or NaN, this sets the `SNESLineSearchReason` to `SNES_LINESEARCH_FAILED_NANORINF` and exits the line search
150 
151 .seealso: [](ch_snes), `SNESSetObjectiveDomainError()`, `PETSC_ERR_CONV_FAILED`, `SNESSetErrorIfNotConverged()`, `SNES_DIVERGED_OBJECTIVE_DOMAIN`, `SNES_DIVERGED_FUNCTION_DOMAIN`,
152           `SNESSetFunctionDomainError()`, `SNESConvergedReason`, `SNES_DIVERGED_OBJECTIVE_NANORINF`, `SNES_DIVERGED_FUNCTION_NAN`, `SNESLineSearchCheckObjectiveDomainError()`
153 MC*/
154 #define SNESLineSearchCheckObjectiveDomainError(snes, fobj) \
155   do { \
156     if (snes->errorifnotconverged) { \
157       PetscCheck(!snes->objectivedomainerror, PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due objective domain error"); \
158       PetscCheck(!PetscIsInfOrNanReal(fobj), PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due to infinity or NaN norm"); \
159     } \
160     if (snes->objectivedomainerror) { \
161       snes->linesearch->reason   = SNES_LINESEARCH_FAILED_OBJECTIVE_DOMAIN; \
162       snes->objectivedomainerror = PETSC_FALSE; \
163       PetscFunctionReturn(PETSC_SUCCESS); \
164     } else if (PetscIsInfOrNanReal(fobj)) { \
165       snes->linesearch->reason = SNES_LINESEARCH_FAILED_NANORINF; \
166       PetscFunctionReturn(PETSC_SUCCESS); \
167     } \
168   } while (0)
169 
170 /*MC
171   SNESLineSearchCheckJacobianDomainError - Called after a `SNESComputeJacobian()` in a SNES line search to check if `SNESSetJacobianDomainError()` has been called.
172 
173   Synopsis:
174   #include <snesimpl.h>
175   void SNESLineSearchCheckJacobian(SNES snes, SNESLineSearch ls)
176 
177   Collective
178 
179   Input Parameters:
180 + snes - the `SNES` solver object
181 - ls   - the `SNESLineSearch` object
182 
183   Level: developer
184 
185   Notes:
186   This turns the non-collective `SNESSetJacobianDomainError()` into a collective operation
187 
188   This check is done in debug mode or if `SNESSetCheckJacobianDomainError()` has been called
189 
190 .seealso: [](ch_snes), `SNESSetCheckJacobianDomainError()`, `SNESSetFunctionDomainError()`, `PETSC_ERR_CONV_FAILED`, `SNESSetErrorIfNotConverged()`, `SNES_DIVERGED_FUNCTION_DOMAIN`,
191           `SNESConvergedReason`, `SNES_DIVERGED_FUNCTION_NAN`
192 MC*/
193 #define SNESLineSearchCheckJacobianDomainError(snes, ls) \
194   do { \
195     if (snes->checkjacdomainerror) { \
196       PetscBool jacobiandomainerror; \
197       PetscCallMPI(MPIU_Allreduce(&snes->jacobiandomainerror, &jacobiandomainerror, 1, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)ls))); \
198       if (jacobiandomainerror) { \
199         ls->reason                = SNES_LINESEARCH_FAILED_JACOBIAN_DOMAIN; \
200         snes->jacobiandomainerror = PETSC_FALSE; \
201         PetscCheck(!snes->errorifnotconverged, PetscObjectComm((PetscObject)ls), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due to Jacobian domain error"); \
202         PetscFunctionReturn(PETSC_SUCCESS); \
203       } \
204     } \
205   } while (0)
206