Actual source code: ls.c

petsc-3.4.5 2014-06-29
  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:   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
137:   PetscReal           fnorm,gnorm,xnorm,ynorm;
138:   Vec                 Y,X,F,G,W,FPC;
139:   KSPConvergedReason  kspreason;
140:   PetscBool           domainerror;
141:   SNESLineSearch      linesearch;
142:   SNESConvergedReason reason;

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

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

156:   PetscObjectAMSTakeAccess((PetscObject)snes);
157:   snes->iter = 0;
158:   snes->norm = 0.0;
159:   PetscObjectAMSGrantAccess((PetscObject)snes);
160:   SNESGetLineSearch(snes, &linesearch);
161:   if (!snes->vec_func_init_set) {
162:     SNESComputeFunction(snes,X,F);
163:     SNESGetFunctionDomainError(snes, &domainerror);
164:     if (domainerror) {
165:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
166:       return(0);
167:     }
168:   } else snes->vec_func_init_set = PETSC_FALSE;

170:   if (!snes->norm_init_set) {
171:     VecNormBegin(F,NORM_2,&fnorm);        /* fnorm <- ||F||  */
172:     VecNormEnd(F,NORM_2,&fnorm);
173:     if (PetscIsInfOrNanReal(fnorm)) {
174:       snes->reason = SNES_DIVERGED_FNORM_NAN;
175:       return(0);
176:     }
177:   } else {
178:     fnorm               = snes->norm_init;
179:     snes->norm_init_set = PETSC_FALSE;
180:   }
181:   PetscObjectAMSTakeAccess((PetscObject)snes);
182:   snes->norm = fnorm;
183:   PetscObjectAMSGrantAccess((PetscObject)snes);
184:   SNESLogConvergenceHistory(snes,fnorm,0);
185:   SNESMonitor(snes,0,fnorm);

187:   /* set parameter for default relative tolerance convergence test */
188:   snes->ttol = fnorm*snes->rtol;
189:   /* test convergence */
190:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
191:   if (snes->reason) return(0);

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

195:     /* Call general purpose update function */
196:     if (snes->ops->update) {
197:       (*snes->ops->update)(snes, snes->iter);
198:     }

200:     /* apply the nonlinear preconditioner if it's right preconditioned */
201:     if (snes->pc && snes->pcside == PC_RIGHT) {
202:       SNESSetInitialFunction(snes->pc, F);
203:       SNESSetInitialFunctionNorm(snes->pc, fnorm);
204:       PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);
205:       SNESSolve(snes->pc, snes->vec_rhs, X);
206:       PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,snes->vec_rhs,0);
207:       SNESGetConvergedReason(snes->pc,&reason);
208:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
209:         snes->reason = SNES_DIVERGED_INNER;
210:         return(0);
211:       }
212:       SNESGetFunction(snes->pc, &FPC, NULL, NULL);
213:       VecCopy(FPC, F);
214:       SNESGetFunctionNorm(snes->pc, &fnorm);
215:     }

217:     /* Solve J Y = F, where J is Jacobian matrix */
218:     SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);
219:     KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);
220:     KSPSolve(snes->ksp,F,Y);
221:     KSPGetConvergedReason(snes->ksp,&kspreason);
222:     if (kspreason < 0) {
223:       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
224:         PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);
225:         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
226:         break;
227:       }
228:     }
229:     KSPGetIterationNumber(snes->ksp,&lits);
230:     snes->linear_its += lits;
231:     PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);

233:     if (PetscLogPrintInfo) {
234:       SNESNEWTONLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);
235:     }

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

287:    Input Parameter:
288: .  snes - the SNES context
289: .  x - the solution vector

291:    Application Interface Routine: SNESSetUp()

293:    Notes:
294:    For basic use of the SNES solvers, the user need not explicitly call
295:    SNESSetUp(), since these actions will automatically occur during
296:    the call to SNESSolve().
297:  */
300: PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
301: {

305:   SNESSetWorkVecs(snes,2);
306:   SNESSetUpMatrices(snes);
307:   return(0);
308: }
309: /* -------------------------------------------------------------------------- */

313: PetscErrorCode SNESReset_NEWTONLS(SNES snes)
314: {
316:   return(0);
317: }

319: /*
320:    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
321:    with SNESCreate_NEWTONLS().

323:    Input Parameter:
324: .  snes - the SNES context

326:    Application Interface Routine: SNESDestroy()
327:  */
330: PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
331: {

335:   SNESReset_NEWTONLS(snes);
336:   PetscFree(snes->data);
337:   return(0);
338: }
339: /* -------------------------------------------------------------------------- */

341: /*
342:    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.

344:    Input Parameters:
345: .  SNES - the SNES context
346: .  viewer - visualization context

348:    Application Interface Routine: SNESView()
349: */
352: static PetscErrorCode SNESView_NEWTONLS(SNES snes,PetscViewer viewer)
353: {
355:   PetscBool      iascii;

358:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
359:   if (iascii) {
360:   }
361:   return(0);
362: }

364: /* -------------------------------------------------------------------------- */
365: /*
366:    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.

368:    Input Parameter:
369: .  snes - the SNES context

371:    Application Interface Routine: SNESSetFromOptions()
372: */
375: static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes)
376: {
378:   SNESLineSearch linesearch;

381:   PetscOptionsHead("SNESNEWTONLS options");
382:   PetscOptionsTail();
383:   /* set the default line search type */
384:   if (!snes->linesearch) {
385:     SNESGetLineSearch(snes, &linesearch);
386:     SNESLineSearchSetType(linesearch, SNESLINESEARCHBT);
387:   }
388:   return(0);
389: }

391: /* -------------------------------------------------------------------------- */
392: /*MC
393:       SNESNEWTONLS - Newton based nonlinear solver that uses a line search

395:    Options Database:
396: +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
397: .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
398: .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch
399: .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
400: .   -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)
401: .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
402: .   -snes_linesearch_monitor - print information about progress of line searches
403: -   -snes_linesearch_damping - damping factor used for basic line search

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

407:    Level: beginner

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

412: M*/
415: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
416: {
418:   SNES_NEWTONLS  *neP;

421:   snes->ops->setup          = SNESSetUp_NEWTONLS;
422:   snes->ops->solve          = SNESSolve_NEWTONLS;
423:   snes->ops->destroy        = SNESDestroy_NEWTONLS;
424:   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
425:   snes->ops->view           = SNESView_NEWTONLS;
426:   snes->ops->reset          = SNESReset_NEWTONLS;

428:   snes->usesksp = PETSC_TRUE;
429:   snes->usespc  = PETSC_TRUE;
430:   PetscNewLog(snes,SNES_NEWTONLS,&neP);
431:   snes->data    = (void*)neP;
432:   return(0);
433: }