xref: /petsc/src/snes/impls/ls/ls.c (revision 0619917b5a674bb687c64e7daba2ab22be99af31)
1 
2 #include <../src/snes/impls/ls/lsimpl.h>
3 
4 /*
5      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
6     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
7     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
8     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
9 */
10 static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin)
11 {
12   PetscReal a1;
13   PetscBool hastranspose;
14   Vec       W;
15   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
16 
17   PetscFunctionBegin;
18   *ismin = PETSC_FALSE;
19   PetscCall(SNESGetObjective(snes, &objective, NULL));
20   if (!objective) {
21     PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
22     PetscCall(VecDuplicate(F, &W));
23     if (hastranspose) {
24       /* Compute || J^T F|| */
25       PetscCall(MatMultTranspose(A, F, W));
26       PetscCall(VecNorm(W, NORM_2, &a1));
27       PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm)));
28       if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
29     } else {
30       Vec         work;
31       PetscScalar result;
32       PetscReal   wnorm;
33 
34       PetscCall(VecSetRandom(W, NULL));
35       PetscCall(VecNorm(W, NORM_2, &wnorm));
36       PetscCall(VecDuplicate(W, &work));
37       PetscCall(MatMult(A, W, work));
38       PetscCall(VecDot(F, work, &result));
39       PetscCall(VecDestroy(&work));
40       a1 = PetscAbsScalar(result) / (fnorm * wnorm);
41       PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1));
42       if (a1 < 1.e-4) *ismin = PETSC_TRUE;
43     }
44     PetscCall(VecDestroy(&W));
45   }
46   PetscFunctionReturn(PETSC_SUCCESS);
47 }
48 
49 /*
50      Checks if J^T(F - J*X) = 0
51 */
52 static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X)
53 {
54   PetscReal a1, a2;
55   PetscBool hastranspose;
56   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
57 
58   PetscFunctionBegin;
59   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
60   PetscCall(SNESGetObjective(snes, &objective, NULL));
61   if (hastranspose && !objective) {
62     Vec W1, W2;
63 
64     PetscCall(VecDuplicate(F, &W1));
65     PetscCall(VecDuplicate(F, &W2));
66     PetscCall(MatMult(A, X, W1));
67     PetscCall(VecAXPY(W1, -1.0, F));
68 
69     /* Compute || J^T W|| */
70     PetscCall(MatMultTranspose(A, W1, W2));
71     PetscCall(VecNorm(W1, NORM_2, &a1));
72     PetscCall(VecNorm(W2, NORM_2, &a2));
73     if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n", (double)(a2 / a1)));
74     PetscCall(VecDestroy(&W1));
75     PetscCall(VecDestroy(&W2));
76   }
77   PetscFunctionReturn(PETSC_SUCCESS);
78 }
79 
80 /*
81      This file implements a truncated Newton method with a line search,
82      for solving a system of nonlinear equations, using the KSP, Vec,
83      and Mat interfaces for linear solvers, vectors, and matrices,
84      respectively.
85 
86      The following basic routines are required for each nonlinear solver:
87           SNESCreate_XXX()          - Creates a nonlinear solver context
88           SNESSetFromOptions_XXX()  - Sets runtime options
89           SNESSolve_XXX()           - Solves the nonlinear system
90           SNESDestroy_XXX()         - Destroys the nonlinear solver context
91      The suffix "_XXX" denotes a particular implementation, in this case
92      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
93      systems of nonlinear equations with a line search (LS) method.
94      These routines are actually called via the common user interface
95      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
96      SNESDestroy(), so the application code interface remains identical
97      for all nonlinear solvers.
98 
99      Another key routine is:
100           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
101      by setting data structures and options.   The interface routine SNESSetUp()
102      is not usually called directly by the user, but instead is called by
103      SNESSolve() if necessary.
104 
105      Additional basic routines are:
106           SNESView_XXX()            - Prints details of runtime options that
107                                       have actually been used.
108      These are called by application codes via the interface routines
109      SNESView().
110 
111      The various types of solvers (preconditioners, Krylov subspace methods,
112      nonlinear solvers, timesteppers) are all organized similarly, so the
113      above description applies to these categories also.
114 
115 */
116 
117 // PetscClangLinter pragma disable: -fdoc-sowing-chars
118 /*
119   SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
120   method with a line search.
121 
122   Input Parameter:
123 . snes - the SNES context
124 
125 */
126 PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
127 {
128   PetscInt             maxits, i, lits;
129   SNESLineSearchReason lssucceed;
130   PetscReal            fnorm, xnorm, ynorm;
131   Vec                  Y, X, F;
132   SNESLineSearch       linesearch;
133   SNESConvergedReason  reason;
134 #if defined(PETSC_USE_INFO)
135   PetscReal gnorm;
136 #endif
137 
138   PetscFunctionBegin;
139   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
140 
141   snes->numFailures            = 0;
142   snes->numLinearSolveFailures = 0;
143   snes->reason                 = SNES_CONVERGED_ITERATING;
144 
145   maxits = snes->max_its;        /* maximum number of iterations */
146   X      = snes->vec_sol;        /* solution vector */
147   F      = snes->vec_func;       /* residual vector */
148   Y      = snes->vec_sol_update; /* newton step */
149 
150   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
151   snes->iter = 0;
152   snes->norm = 0.0;
153   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
154   PetscCall(SNESGetLineSearch(snes, &linesearch));
155 
156   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
157   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
158     PetscCall(SNESApplyNPC(snes, X, NULL, F));
159     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
160     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
161       PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
162       PetscFunctionReturn(PETSC_SUCCESS);
163     }
164 
165     PetscCall(VecNormBegin(F, NORM_2, &fnorm));
166     PetscCall(VecNormEnd(F, NORM_2, &fnorm));
167   } else {
168     if (!snes->vec_func_init_set) {
169       PetscCall(SNESComputeFunction(snes, X, F));
170     } else snes->vec_func_init_set = PETSC_FALSE;
171   }
172 
173   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
174   SNESCheckFunctionNorm(snes, fnorm);
175   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
176   snes->norm = fnorm;
177   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
178   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
179 
180   /* test convergence */
181   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
182   PetscCall(SNESMonitor(snes, 0, fnorm));
183   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
184 
185   for (i = 0; i < maxits; i++) {
186     /* Call general purpose update function */
187     PetscTryTypeMethod(snes, update, snes->iter);
188 
189     /* apply the nonlinear preconditioner */
190     if (snes->npc) {
191       if (snes->npcside == PC_RIGHT) {
192         PetscCall(SNESSetInitialFunction(snes->npc, F));
193         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
194         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
195         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
196         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
197         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
198           PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
199           PetscFunctionReturn(PETSC_SUCCESS);
200         }
201         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
202       } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
203         PetscCall(SNESApplyNPC(snes, X, F, F));
204         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
205         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
206           PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
207           PetscFunctionReturn(PETSC_SUCCESS);
208         }
209       }
210     }
211 
212     /* Solve J Y = F, where J is Jacobian matrix */
213     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
214     SNESCheckJacobianDomainerror(snes);
215     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
216     PetscCall(KSPSolve(snes->ksp, F, Y));
217     SNESCheckKSPSolve(snes);
218     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
219     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
220 
221     if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y));
222 
223 #if defined(PETSC_USE_INFO)
224     gnorm = fnorm;
225 #endif
226     /* Compute a (scaled) negative update in the line search routine:
227          X <- X - lambda*Y
228        and evaluate F = function(X) (depends on the line search).
229     */
230     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y));
231     PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed));
232     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
233     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed));
234     if (snes->reason) break;
235     SNESCheckFunctionNorm(snes, fnorm);
236     if (lssucceed) {
237       if (snes->stol * xnorm > ynorm) {
238         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
239         PetscFunctionReturn(PETSC_SUCCESS);
240       }
241       if (++snes->numFailures >= snes->maxFailures) {
242         PetscBool ismin;
243         snes->reason = SNES_DIVERGED_LINE_SEARCH;
244         PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin));
245         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
246         if (snes->errorifnotconverged && snes->reason) {
247           PetscViewer monitor;
248           PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
249           PetscCheck(monitor, PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due to %s. Suggest running with -snes_linesearch_monitor", SNESConvergedReasons[snes->reason]);
250           SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]);
251         }
252         break;
253       }
254     }
255     /* Monitor convergence */
256     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
257     snes->iter  = i + 1;
258     snes->norm  = fnorm;
259     snes->ynorm = ynorm;
260     snes->xnorm = xnorm;
261     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
262     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
263     /* Test for convergence */
264     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
265     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
266     if (snes->reason) break;
267   }
268   PetscFunctionReturn(PETSC_SUCCESS);
269 }
270 
271 /*
272    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
273    of the SNESNEWTONLS nonlinear solver.
274 
275    Input Parameter:
276 .  snes - the SNES context
277 .  x - the solution vector
278 
279    Application Interface Routine: SNESSetUp()
280 
281  */
282 PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
283 {
284   PetscFunctionBegin;
285   PetscCall(SNESSetUpMatrices(snes));
286   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
287   PetscFunctionReturn(PETSC_SUCCESS);
288 }
289 
290 PetscErrorCode SNESReset_NEWTONLS(SNES snes)
291 {
292   PetscFunctionBegin;
293   PetscFunctionReturn(PETSC_SUCCESS);
294 }
295 
296 /*
297    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
298    with SNESCreate_NEWTONLS().
299 
300    Input Parameter:
301 .  snes - the SNES context
302 
303    Application Interface Routine: SNESDestroy()
304  */
305 PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
306 {
307   PetscFunctionBegin;
308   PetscCall(SNESReset_NEWTONLS(snes));
309   PetscCall(PetscFree(snes->data));
310   PetscFunctionReturn(PETSC_SUCCESS);
311 }
312 
313 /*
314    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
315 
316    Input Parameters:
317 .  SNES - the SNES context
318 .  viewer - visualization context
319 
320    Application Interface Routine: SNESView()
321 */
322 static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer)
323 {
324   PetscBool iascii;
325 
326   PetscFunctionBegin;
327   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
328   if (iascii) { }
329   PetscFunctionReturn(PETSC_SUCCESS);
330 }
331 
332 /*
333    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
334 
335    Input Parameter:
336 .  snes - the SNES context
337 
338    Application Interface Routine: SNESSetFromOptions()
339 */
340 static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject)
341 {
342   PetscFunctionBegin;
343   PetscFunctionReturn(PETSC_SUCCESS);
344 }
345 
346 /*MC
347       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
348 
349    Options Database Keys:
350 +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
351 .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
352 .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`)
353 .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
354 .   -snes_linesearch_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
355 .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
356 .   -snes_linesearch_monitor - print information about progress of line searches
357 -   -snes_linesearch_damping - damping factor used for basic line search
358 
359     Note:
360     This is the default nonlinear solver in `SNES`
361 
362    Level: beginner
363 
364 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()`
365           `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()`
366 M*/
367 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
368 {
369   SNES_NEWTONLS *neP;
370   SNESLineSearch linesearch;
371 
372   PetscFunctionBegin;
373   snes->ops->setup          = SNESSetUp_NEWTONLS;
374   snes->ops->solve          = SNESSolve_NEWTONLS;
375   snes->ops->destroy        = SNESDestroy_NEWTONLS;
376   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
377   snes->ops->view           = SNESView_NEWTONLS;
378   snes->ops->reset          = SNESReset_NEWTONLS;
379 
380   snes->npcside = PC_RIGHT;
381   snes->usesksp = PETSC_TRUE;
382   snes->usesnpc = PETSC_TRUE;
383 
384   PetscCall(SNESGetLineSearch(snes, &linesearch));
385   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
386 
387   snes->alwayscomputesfinalresidual = PETSC_TRUE;
388 
389   PetscCall(PetscNew(&neP));
390   snes->data = (void *)neP;
391   PetscFunctionReturn(PETSC_SUCCESS);
392 }
393