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