Actual source code: ls.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <../src/snes/impls/ls/lsimpl.h>

  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: */
 12: PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
 13: {
 14:   PetscReal      a1;
 16:   PetscBool      hastranspose;

 19:   *ismin = PETSC_FALSE;
 20:   MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
 21:   if (hastranspose) {
 22:     /* Compute || J^T F|| */
 23:     MatMultTranspose(A,F,W);
 24:     VecNorm(W,NORM_2,&a1);
 25:     PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));
 26:     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
 27:   } else {
 28:     Vec         work;
 29:     PetscScalar result;
 30:     PetscReal   wnorm;

 32:     VecSetRandom(W,NULL);
 33:     VecNorm(W,NORM_2,&wnorm);
 34:     VecDuplicate(W,&work);
 35:     MatMult(A,W,work);
 36:     VecDot(F,work,&result);
 37:     VecDestroy(&work);
 38:     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
 39:     PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);
 40:     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
 41:   }
 42:   return(0);
 43: }

 45: /*
 46:      Checks if J^T(F - J*X) = 0
 47: */
 50: PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
 51: {
 52:   PetscReal      a1,a2;
 54:   PetscBool      hastranspose;

 57:   MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);
 58:   if (hastranspose) {
 59:     MatMult(A,X,W1);
 60:     VecAXPY(W1,-1.0,F);

 62:     /* Compute || J^T W|| */
 63:     MatMultTranspose(A,W1,W2);
 64:     VecNorm(W1,NORM_2,&a1);
 65:     VecNorm(W2,NORM_2,&a2);
 66:     if (a1 != 0.0) {
 67:       PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));
 68:     }
 69:   }
 70:   return(0);
 71: }

 73: /*  --------------------------------------------------------------------

 75:      This file implements a truncated Newton method with a line search,
 76:      for solving a system of nonlinear equations, using the KSP, Vec,
 77:      and Mat interfaces for linear solvers, vectors, and matrices,
 78:      respectively.

 80:      The following basic routines are required for each nonlinear solver:
 81:           SNESCreate_XXX()          - Creates a nonlinear solver context
 82:           SNESSetFromOptions_XXX()  - Sets runtime options
 83:           SNESSolve_XXX()           - Solves the nonlinear system
 84:           SNESDestroy_XXX()         - Destroys the nonlinear solver context
 85:      The suffix "_XXX" denotes a particular implementation, in this case
 86:      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
 87:      systems of nonlinear equations with a line search (LS) method.
 88:      These routines are actually called via the common user interface
 89:      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
 90:      SNESDestroy(), so the application code interface remains identical
 91:      for all nonlinear solvers.

 93:      Another key routine is:
 94:           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
 95:      by setting data structures and options.   The interface routine SNESSetUp()
 96:      is not usually called directly by the user, but instead is called by
 97:      SNESSolve() if necessary.

 99:      Additional basic routines are:
100:           SNESView_XXX()            - Prints details of runtime options that
101:                                       have actually been used.
102:      These are called by application codes via the interface routines
103:      SNESView().

105:      The various types of solvers (preconditioners, Krylov subspace methods,
106:      nonlinear solvers, timesteppers) are all organized similarly, so the
107:      above description applies to these categories also.

109:     -------------------------------------------------------------------- */
110: /*
111:    SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
112:    method with a line search.

114:    Input Parameters:
115: .  snes - the SNES context

117:    Output Parameter:
118: .  outits - number of iterations until termination

120:    Application Interface Routine: SNESSolve()

122:    Notes:
123:    This implements essentially a truncated Newton method with a
124:    line search.  By default a cubic backtracking line search
125:    is employed, as described in the text "Numerical Methods for
126:    Unconstrained Optimization and Nonlinear Equations" by Dennis
127:    and Schnabel.
128: */
131: PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
132: {
133:   PetscErrorCode      ierr;
134:   PetscInt            maxits,i,lits;
135:   PetscBool           lssucceed;
136:   PetscReal           fnorm,gnorm,xnorm,ynorm;
137:   Vec                 Y,X,F,G,W;
138:   KSPConvergedReason  kspreason;
139:   PetscBool           domainerror;
140:   SNESLineSearch      linesearch;
141:   SNESConvergedReason reason;

144:   snes->numFailures            = 0;
145:   snes->numLinearSolveFailures = 0;
146:   snes->reason                 = SNES_CONVERGED_ITERATING;

148:   maxits = snes->max_its;               /* maximum number of iterations */
149:   X      = snes->vec_sol;               /* solution vector */
150:   F      = snes->vec_func;              /* residual vector */
151:   Y      = snes->vec_sol_update;        /* newton step */
152:   G      = snes->work[0];
153:   W      = snes->work[1];

155:   PetscObjectSAWsTakeAccess((PetscObject)snes);
156:   snes->iter = 0;
157:   snes->norm = 0.0;
158:   PetscObjectSAWsGrantAccess((PetscObject)snes);
159:   SNESGetLineSearch(snes, &linesearch);

161:   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
162:   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
163:     SNESApplyNPC(snes,X,NULL,F);
164:     SNESGetConvergedReason(snes->pc,&reason);
165:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
166:       snes->reason = SNES_DIVERGED_INNER;
167:       return(0);
168:     }

170:     VecNormBegin(F,NORM_2,&fnorm);
171:     VecNormEnd(F,NORM_2,&fnorm);
172:   } else {
173:     if (!snes->vec_func_init_set) {
174:       SNESComputeFunction(snes,X,F);
175:       SNESGetFunctionDomainError(snes, &domainerror);
176:       if (domainerror) {
177:         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
178:         return(0);
179:       }
180:     } else snes->vec_func_init_set = PETSC_FALSE;
181:   }

183:   VecNorm(F,NORM_2,&fnorm);        /* fnorm <- ||F||  */
184:   if (PetscIsInfOrNanReal(fnorm)) {
185:     snes->reason = SNES_DIVERGED_FNORM_NAN;
186:     return(0);
187:   }

189:   PetscObjectSAWsTakeAccess((PetscObject)snes);
190:   snes->norm = fnorm;
191:   PetscObjectSAWsGrantAccess((PetscObject)snes);
192:   SNESLogConvergenceHistory(snes,fnorm,0);
193:   SNESMonitor(snes,0,fnorm);

195:   /* test convergence */
196:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
197:   if (snes->reason) return(0);

199:   for (i=0; i<maxits; i++) {

201:     /* Call general purpose update function */
202:     if (snes->ops->update) {
203:       (*snes->ops->update)(snes, snes->iter);
204:     }

206:     /* apply the nonlinear preconditioner */
207:     if (snes->pc) {
208:       if (snes->pcside == PC_RIGHT) {
209:         SNESSetInitialFunction(snes->pc, F);
210:         PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);
211:         SNESSolve(snes->pc, snes->vec_rhs, X);
212:         PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);
213:         SNESGetConvergedReason(snes->pc,&reason);
214:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
215:           snes->reason = SNES_DIVERGED_INNER;
216:           return(0);
217:         }
218:         SNESGetNPCFunction(snes,F,&fnorm);
219:       } else if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
220:         SNESApplyNPC(snes,X,F,F);
221:         SNESGetConvergedReason(snes->pc,&reason);
222:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
223:           snes->reason = SNES_DIVERGED_INNER;
224:           return(0);
225:         }
226:       }
227:     }

229:     /* Solve J Y = F, where J is Jacobian matrix */
230:     SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
231:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
232:     KSPSolve(snes->ksp,F,Y);
233:     KSPGetConvergedReason(snes->ksp,&kspreason);
234:     if (kspreason < 0) {
235:       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
236:         PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
237:         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
238:         break;
239:       }
240:     }
241:     KSPGetIterationNumber(snes->ksp,&lits);
242:     snes->linear_its += lits;
243:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);

245:     if (PetscLogPrintInfo) {
246:       SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
247:     }

249:     /* Compute a (scaled) negative update in the line search routine:
250:          X <- X - lambda*Y
251:        and evaluate F = function(X) (depends on the line search).
252:     */
253:     gnorm = fnorm;
254:     SNESLineSearchApply(linesearch, X, F, &fnorm, Y);
255:     SNESLineSearchGetSuccess(linesearch, &lssucceed);
256:     SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
257:     PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)gnorm,(double)fnorm,(double)ynorm,(int)lssucceed);
258:     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
259:     SNESGetFunctionDomainError(snes, &domainerror);
260:     if (domainerror) {
261:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
262:       return(0);
263:     }
264:     if (!lssucceed) {
265:       if (snes->stol*xnorm > ynorm) {
266:         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
267:         return(0);
268:       }
269:       if (++snes->numFailures >= snes->maxFailures) {
270:         PetscBool ismin;
271:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
272:         SNESNEWTONLSCheckLocalMin_Private(snes,snes->jacobian,F,W,fnorm,&ismin);
273:         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
274:         break;
275:       }
276:     }
277:     /* Monitor convergence */
278:     PetscObjectSAWsTakeAccess((PetscObject)snes);
279:     snes->iter = i+1;
280:     snes->norm = fnorm;
281:     PetscObjectSAWsGrantAccess((PetscObject)snes);
282:     SNESLogConvergenceHistory(snes,snes->norm,lits);
283:     SNESMonitor(snes,snes->iter,snes->norm);
284:     /* Test for convergence */
285:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
286:     if (snes->reason) break;
287:   }
288:   if (i == maxits) {
289:     PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);
290:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
291:   }
292:   return(0);
293: }
294: /* -------------------------------------------------------------------------- */
295: /*
296:    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
297:    of the SNESNEWTONLS nonlinear solver.

299:    Input Parameter:
300: .  snes - the SNES context
301: .  x - the solution vector

303:    Application Interface Routine: SNESSetUp()

305:    Notes:
306:    For basic use of the SNES solvers, the user need not explicitly call
307:    SNESSetUp(), since these actions will automatically occur during
308:    the call to SNESSolve().
309:  */
312: PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
313: {

317:   SNESSetWorkVecs(snes,2);
318:   SNESSetUpMatrices(snes);
319:   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
320:   return(0);
321: }
322: /* -------------------------------------------------------------------------- */

326: PetscErrorCode SNESReset_NEWTONLS(SNES snes)
327: {
329:   return(0);
330: }

332: /*
333:    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
334:    with SNESCreate_NEWTONLS().

336:    Input Parameter:
337: .  snes - the SNES context

339:    Application Interface Routine: SNESDestroy()
340:  */
343: PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
344: {

348:   SNESReset_NEWTONLS(snes);
349:   PetscFree(snes->data);
350:   return(0);
351: }
352: /* -------------------------------------------------------------------------- */

354: /*
355:    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.

357:    Input Parameters:
358: .  SNES - the SNES context
359: .  viewer - visualization context

361:    Application Interface Routine: SNESView()
362: */
365: static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer)
366: {
368:   PetscBool      iascii;

371:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
372:   if (iascii) {
373:   }
374:   return(0);
375: }

377: /* -------------------------------------------------------------------------- */
378: /*
379:    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.

381:    Input Parameter:
382: .  snes - the SNES context

384:    Application Interface Routine: SNESSetFromOptions()
385: */
388: static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes)
389: {
391:   SNESLineSearch linesearch;

394:   if (!snes->linesearch) {
395:     SNESGetLineSearch(snes, &linesearch);
396:     SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);
397:   }
398:   return(0);
399: }

401: /* -------------------------------------------------------------------------- */
402: /*MC
403:       SNESNEWTONLS - Newton based nonlinear solver that uses a line search

405:    Options Database:
406: +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
407: .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
408: .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
409: .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
410: .   -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)
411: .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
412: .   -snes_linesearch_monitor - print information about progress of line searches
413: -   -snes_linesearch_damping - damping factor used for basic line search

415:     Notes: This is the default nonlinear solver in SNES

417:    Level: beginner

419: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONTR, SNESQN, SNESLineSearchSetType(), SNESLineSearchSetOrder()
420:            SNESLineSearchSetPostCheck(), SNESLineSearchSetPreCheck() SNESLineSearchSetComputeNorms()

422: M*/
425: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
426: {
428:   SNES_NEWTONLS  *neP;

431:   snes->ops->setup          = SNESSetUp_NEWTONLS;
432:   snes->ops->solve          = SNESSolve_NEWTONLS;
433:   snes->ops->destroy        = SNESDestroy_NEWTONLS;
434:   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
435:   snes->ops->view           = SNESView_NEWTONLS;
436:   snes->ops->reset          = SNESReset_NEWTONLS;

438:   snes->pcside  = PC_RIGHT;
439:   snes->usesksp = PETSC_TRUE;
440:   snes->usespc  = PETSC_TRUE;
441:   PetscNewLog(snes,&neP);
442:   snes->data    = (void*)neP;
443:   return(0);
444: }