Actual source code: iterativ.c

  1: #define PETSCKSP_DLL

  3: /*
  4:    This file contains some simple default routines.  
  5:    These routines should be SHORT, since they will be included in every
  6:    executable image that uses the iterative routines (note that, through
  7:    the registry system, we provide a way to load only the truely necessary
  8:    files) 
  9:  */
 10:  #include private/kspimpl.h

 14: /*
 15:   KSPDefaultFreeWork - Free work vectors

 17:   Input Parameters:
 18: . ksp  - iterative context
 19:  */
 20: PetscErrorCode KSPDefaultFreeWork(KSP ksp)
 21: {
 25:   if (ksp->work)  {
 26:     VecDestroyVecs(ksp->work,ksp->nwork);
 27:     ksp->work = PETSC_NULL;
 28:   }
 29:   return(0);
 30: }

 34: /*@
 35:    KSPGetResidualNorm - Gets the last (approximate preconditioned)
 36:    residual norm that has been computed.
 37:  
 38:    Not Collective

 40:    Input Parameters:
 41: .  ksp - the iterative context

 43:    Output Parameters:
 44: .  rnorm - residual norm

 46:    Level: intermediate

 48: .keywords: KSP, get, residual norm

 50: .seealso: KSPBuildResidual()
 51: @*/
 52: PetscErrorCode  KSPGetResidualNorm(KSP ksp,PetscReal *rnorm)
 53: {
 57:   *rnorm = ksp->rnorm;
 58:   return(0);
 59: }

 63: /*@
 64:    KSPGetIterationNumber - Gets the current iteration number; if the 
 65:          KSPSolve() is complete, returns the number of iterations
 66:          used.
 67:  
 68:    Not Collective

 70:    Input Parameters:
 71: .  ksp - the iterative context

 73:    Output Parameters:
 74: .  its - number of iterations

 76:    Level: intermediate

 78:    Notes:
 79:       During the ith iteration this returns i-1
 80: .keywords: KSP, get, residual norm

 82: .seealso: KSPBuildResidual(), KSPGetResidualNorm()
 83: @*/
 84: PetscErrorCode  KSPGetIterationNumber(KSP ksp,PetscInt *its)
 85: {
 89:   *its = ksp->its;
 90:   return(0);
 91: }

 95: /*@C
 96:     KSPMonitorSingularValue - Prints the two norm of the true residual and
 97:     estimation of the extreme singular values of the preconditioned problem
 98:     at each iteration.
 99:  
100:     Collective on KSP

102:     Input Parameters:
103: +   ksp - the iterative context
104: .   n  - the iteration
105: -   rnorm - the two norm of the residual

107:     Options Database Key:
108: .   -ksp_monitor_singular_value - Activates KSPMonitorSingularValue()

110:     Notes:
111:     The CG solver uses the Lanczos technique for eigenvalue computation, 
112:     while GMRES uses the Arnoldi technique; other iterative methods do
113:     not currently compute singular values.

115:     Level: intermediate

117: .keywords: KSP, CG, default, monitor, extreme, singular values, Lanczos, Arnoldi

119: .seealso: KSPComputeExtremeSingularValues()
120: @*/
121: PetscErrorCode  KSPMonitorSingularValue(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
122: {
123:   PetscReal               emin,emax,c;
124:   PetscErrorCode          ierr;
125:   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;

129:   if (!dummy) {PetscViewerASCIIMonitorCreate(((PetscObject)ksp)->comm,"stdout",0,&viewer);}
130:   if (!ksp->calc_sings) {
131:     PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,rnorm);
132:   } else {
133:     KSPComputeExtremeSingularValues(ksp,&emax,&emin);
134:     c = emax/emin;
135:     PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP Residual norm %14.12e %% max %G min %G max/min %G\n",n,rnorm,emax,emin,c);
136:   }
137:   if (!dummy) {PetscViewerASCIIMonitorDestroy(viewer);}
138:   return(0);
139: }

143: /*@C
144:    KSPMonitorSolution - Monitors progress of the KSP solvers by calling 
145:    VecView() for the approximate solution at each iteration.

147:    Collective on KSP

149:    Input Parameters:
150: +  ksp - the KSP context
151: .  its - iteration number
152: .  fgnorm - 2-norm of residual (or gradient)
153: -  dummy - either a viewer or PETSC_NULL

155:    Level: intermediate

157:    Notes:
158:     For some Krylov methods such as GMRES constructing the solution at
159:   each iteration is expensive, hence using this will slow the code.

161: .keywords: KSP, nonlinear, vector, monitor, view

163: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView()
164: @*/
165: PetscErrorCode  KSPMonitorSolution(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
166: {
168:   Vec            x;
169:   PetscViewer    viewer = (PetscViewer) dummy;

172:   KSPBuildSolution(ksp,PETSC_NULL,&x);
173:   if (!viewer) {
174:     MPI_Comm comm;
175:     PetscObjectGetComm((PetscObject)ksp,&comm);
176:     viewer = PETSC_VIEWER_DRAW_(comm);
177:   }
178:   VecView(x,viewer);

180:   return(0);
181: }

185: /*@C
186:    KSPMonitorDefault - Print the residual norm at each iteration of an
187:    iterative solver.

189:    Collective on KSP

191:    Input Parameters:
192: +  ksp   - iterative context
193: .  n     - iteration number
194: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
195: -  dummy - unused monitor context 

197:    Level: intermediate

199: .keywords: KSP, default, monitor, residual

201: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGCreate()
202: @*/
203: PetscErrorCode  KSPMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
204: {
205:   PetscErrorCode          ierr;
206:   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;

209:   if (!dummy) {PetscViewerASCIIMonitorCreate(((PetscObject)ksp)->comm,"stdout",0,&viewer);}
210:   if (n == 0 && ((PetscObject)ksp)->prefix) {
211:     PetscViewerASCIIMonitorPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
212:   }
213:   PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,rnorm);
214:   if (!dummy) {PetscViewerASCIIMonitorDestroy(viewer);}
215:   return(0);
216: }

220: /*@C
221:    KSPMonitorTrueResidualNorm - Prints the true residual norm as well as the preconditioned
222:    residual norm at each iteration of an iterative solver.

224:    Collective on KSP

226:    Input Parameters:
227: +  ksp   - iterative context
228: .  n     - iteration number
229: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
230: -  dummy - unused monitor context 

232:    Options Database Key:
233: .  -ksp_monitor_true_residual - Activates KSPMonitorTrueResidualNorm()

235:    Notes:
236:    When using right preconditioning, these values are equivalent.

238:    When using either ICC or ILU preconditioners in BlockSolve95 
239:    (via MATMPIROWBS matrix format), then use this monitor will
240:    print both the residual norm associated with the original
241:    (unscaled) matrix.

243:    Level: intermediate

245: .keywords: KSP, default, monitor, residual

247: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGCreate()
248: @*/
249: PetscErrorCode  KSPMonitorTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
250: {
251:   PetscErrorCode          ierr;
252:   Vec                     resid,work;
253:   PetscReal               scnorm,bnorm;
254:   PC                      pc;
255:   Mat                     A,B;
256:   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;
257: 
259:   if (!dummy) {PetscViewerASCIIMonitorCreate(((PetscObject)ksp)->comm,"stdout",0,&viewer);}
260:   VecDuplicate(ksp->vec_rhs,&work);
261:   KSPBuildResidual(ksp,0,work,&resid);

263:   /*
264:      Unscale the residual if the matrix is, for example, a BlockSolve matrix
265:     but only if both matrices are the same matrix, since only then would 
266:     they be scaled.
267:   */
268:   VecCopy(resid,work);
269:   KSPGetPC(ksp,&pc);
270:   PCGetOperators(pc,&A,&B,PETSC_NULL);
271:   if (A == B) {
272:     MatUnScaleSystem(A,work,PETSC_NULL);
273:   }
274:   VecNorm(work,NORM_2,&scnorm);
275:   VecDestroy(work);
276:   VecNorm(ksp->vec_rhs,NORM_2,&bnorm);
277:   PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP preconditioned resid norm %14.12e true resid norm %14.12e ||Ae||/||Ax|| %14.12e\n",n,rnorm,scnorm,scnorm/bnorm);
278:   if (!dummy) {PetscViewerASCIIMonitorDestroy(viewer);}
279:   return(0);
280: }

284: PetscErrorCode  KSPMonitorRange_Private(KSP ksp,PetscInt it,PetscReal *per)
285: {
286:   PetscErrorCode          ierr;
287:   Vec                     resid,work;
288:   PetscReal               rmax,pwork;
289:   PC                      pc;
290:   Mat                     A,B;
291:   PetscInt                i,n,N;
292:   PetscScalar             *r;

295:   VecDuplicate(ksp->vec_rhs,&work);
296:   KSPBuildResidual(ksp,0,work,&resid);

298:   /*
299:      Unscale the residual if the matrix is, for example, a BlockSolve matrix
300:     but only if both matrices are the same matrix, since only then would 
301:     they be scaled.
302:   */
303:   VecCopy(resid,work);
304:   KSPGetPC(ksp,&pc);
305:   PCGetOperators(pc,&A,&B,PETSC_NULL);
306:   if (A == B) {
307:     MatUnScaleSystem(A,work,PETSC_NULL);
308:   }
309:   VecNorm(work,NORM_INFINITY,&rmax);
310:   VecGetLocalSize(work,&n);
311:   VecGetSize(work,&N);
312:   VecGetArray(work,&r);
313:   pwork = 0.0;
314:   for (i=0; i<n; i++) {
315:     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
316:   }
317:   MPI_Allreduce(&pwork,per,1,MPIU_REAL,PetscSum_Op,((PetscObject)ksp)->comm);
318:   VecRestoreArray(work,&r);
319:   VecDestroy(work);
320:   *per  = *per/N;
321:   return(0);
322: }

326: /*@C
327:    KSPMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.

329:    Collective on KSP

331:    Input Parameters:
332: +  ksp   - iterative context
333: .  it    - iteration number
334: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
335: -  dummy - unused monitor context 

337:    Options Database Key:
338: .  -ksp_monitor_range - Activates KSPMonitorRange()


341:    Notes:
342:    When using either ICC or ILU preconditioners in BlockSolve95 
343:    (via MATMPIROWBS matrix format), then use this monitor will
344:    print both the residual norm associated with the original
345:    (unscaled) matrix.

347:    Level: intermediate

349: .keywords: KSP, default, monitor, residual

351: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGCreate()
352: @*/
353: PetscErrorCode  KSPMonitorRange(KSP ksp,PetscInt it,PetscReal rnorm,void *dummy)
354: {
355:   PetscErrorCode          ierr;
356:   PetscReal               perc,rel;
357:   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;
358:   /* should be in a MonitorRangeContext */
359:   static PetscReal        prev;

362:   if (!it) prev = rnorm;
363:   if (!dummy) {PetscViewerASCIIMonitorCreate(((PetscObject)ksp)->comm,"stdout",0,&viewer);}
364:   KSPMonitorRange_Private(ksp,it,&perc);

366:   rel  = (prev - rnorm)/prev;
367:   prev = rnorm;
368:   PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e \n",it,rnorm,100.0*perc,rel,rel/perc);
369:   if (!dummy) {PetscViewerASCIIMonitorDestroy(viewer);}
370:   return(0);
371: }

375: /*
376:   Default (short) KSP Monitor, same as KSPMonitorDefault() except
377:   it prints fewer digits of the residual as the residual gets smaller.
378:   This is because the later digits are meaningless and are often 
379:   different on different machines; by using this routine different 
380:   machines will usually generate the same output.
381: */
382: PetscErrorCode  KSPMonitorDefaultShort(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy)
383: {
384:   PetscErrorCode          ierr;
385:   PetscViewerASCIIMonitor viewer = (PetscViewerASCIIMonitor) dummy;

388:   if (!dummy) {PetscViewerASCIIMonitorCreate(((PetscObject)ksp)->comm,"stdout",0,&viewer);}

390:   if (fnorm > 1.e-9) {
391:     PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP Residual norm %G \n",its,fnorm);
392:   } else if (fnorm > 1.e-11){
393:     PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP Residual norm %5.3e \n",its,fnorm);
394:   } else {
395:     PetscViewerASCIIMonitorPrintf(viewer,"%3D KSP Residual norm < 1.e-11\n",its);
396:   }
397:   if (!dummy) {PetscViewerASCIIMonitorDestroy(viewer);}
398:   return(0);
399: }

403: /*@C
404:    KSPSkipConverged - Convergence test that do not return as converged
405:    until the maximum number of iterations is reached.

407:    Collective on KSP

409:    Input Parameters:
410: +  ksp   - iterative context
411: .  n     - iteration number
412: .  rnorm - 2-norm residual value (may be estimated)
413: -  dummy - unused convergence context 

415:    Returns:
416: .  reason - KSP_CONVERGED_ITERATING, KSP_CONVERGED_ITS

418:    Notes: 
419:    This should be used as the convergence test with the option
420:    KSPSetNormType(ksp,KSP_NORM_NO), since norms of the residual are
421:    not computed. Convergence is then declared after the maximum number
422:    of iterations have been reached. Useful when one is using CG or
423:    BiCGStab as a smoother.
424:                     
425:    Level: advanced

427: .keywords: KSP, default, convergence, residual

429: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSetNormType()
430: @*/
431: PetscErrorCode  KSPSkipConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
432: {
436:   *reason = KSP_CONVERGED_ITERATING;
437:   if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
438:   return(0);
439: }

441: typedef struct {
442:   PetscTruth initialrtol;    /* default relative residual decrease is computing from initial residual, not rhs */
443:   PetscTruth mininitialrtol; /* default relative residual decrease is computing from min of initial residual and rhs */
444:   Vec        work;
445: } KSPDefaultConvergedCtx;

449: /*@C
450:    KSPDefaultConvergedCreate - Creates and initializes the space used by the KSPDefaultConverged() function context

452:    Collective on KSP

454:    Output Parameter:
455: .  ctx - convergence context 

457:    Level: intermediate

459: .keywords: KSP, default, convergence, residual

461: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(),
462:           KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm(), KSPDefaultConvergedDestroy()
463: @*/
464: PetscErrorCode  KSPDefaultConvergedCreate(void **ctx)
465: {
466:   PetscErrorCode         ierr;
467:   KSPDefaultConvergedCtx *cctx;

470:   PetscNew(KSPDefaultConvergedCtx,&cctx);
471:   *ctx = cctx;
472:   return(0);
473: }

477: /*@
478:    KSPDefaultConvergedSetUIRNorm - makes the default convergence test use || B*(b - A*(initial guess))||
479:       instead of || B*b ||. In the case of right preconditioner or if KSPSetNormType(ksp,KSP_NORM_UNPRECONDIITONED)
480:       is used there is no B in the above formula. UIRNorm is short for Use Initial Residual Norm.

482:    Collective on KSP

484:    Input Parameters:
485: .  ksp   - iterative context

487:    Options Database:
488: .   -ksp_converged_use_initial_residual_norm

490:    Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.

492:    The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
493:    are defined in petscksp.h.

495:    If the convergence test is not KSPDefaultConverged() then this is ignored.
496:    Level: intermediate

498: .keywords: KSP, default, convergence, residual

500: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUMIRNorm()
501: @*/
502: PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP ksp)
503: {
504:   KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) ksp->cnvP;

508:   if (ksp->converged != KSPDefaultConverged) return(0);
509:   if (ctx->mininitialrtol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPDefaultConvergedSetUIRNorm() and KSPDefaultConvergedSetUMIRNorm() together");
510:   ctx->initialrtol = PETSC_TRUE;
511:   return(0);
512: }

516: /*@
517:    KSPDefaultConvergedSetUMIRNorm - makes the default convergence test use min(|| B*(b - A*(initial guess))||,|| B*b ||)
518:       In the case of right preconditioner or if KSPSetNormType(ksp,KSP_NORM_UNPRECONDIITONED)
519:       is used there is no B in the above formula. UMIRNorm is short for Use Minimum Initial Residual Norm.

521:    Collective on KSP

523:    Input Parameters:
524: .  ksp   - iterative context

526:    Options Database:
527: .   -ksp_converged_use_min_initial_residual_norm

529:    Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.

531:    The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
532:    are defined in petscksp.h.

534:    Level: intermediate

536: .keywords: KSP, default, convergence, residual

538: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm()
539: @*/
540: PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP ksp)
541: {
542:   KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) ksp->cnvP;

546:   if (ksp->converged != KSPDefaultConverged) return(0);
547:   if (ctx->initialrtol) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPDefaultConvergedSetUIRNorm() and KSPDefaultConvergedSetUMIRNorm() together");
548:   ctx->mininitialrtol = PETSC_TRUE;
549:    return(0);
550: }

554: /*@C
555:    KSPDefaultConverged - Determines convergence of
556:    the iterative solvers (default code).

558:    Collective on KSP

560:    Input Parameters:
561: +  ksp   - iterative context
562: .  n     - iteration number
563: .  rnorm - 2-norm residual value (may be estimated)
564: -  dummy - unused convergence context 

566:    reason is set to:
567: +   positive - if the iteration has converged;
568: .   negative - if residual norm exceeds divergence threshold;
569: -   0 - otherwise.

571:    Notes:
572:    KSPDefaultConverged() reaches convergence when
573: $      rnorm < MAX (rtol * rnorm_0, abstol);
574:    Divergence is detected if
575: $      rnorm > dtol * rnorm_0,

577:    where 
578: +     rtol = relative tolerance,
579: .     abstol = absolute tolerance.
580: .     dtol = divergence tolerance,
581: -     rnorm_0 is the two norm of the right hand side. When initial guess is non-zero you 
582:           can call KSPDefaultConvergedSetUIRNorm() to use the norm of (b - A*(initial guess))
583:           as the starting point for relative norm convergence testing.

585:    Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.

587:    The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
588:    are defined in petscksp.h.

590:    Level: intermediate

592: .keywords: KSP, default, convergence, residual

594: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(),
595:           KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm(), KSPDefaultConvergedCreate(), KSPDefaultConvergedDestroy()
596: @*/
597: PetscErrorCode  KSPDefaultConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
598: {
599:   PetscErrorCode         ierr;
600:   KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) dummy;
601:   KSPNormType            normtype;

606:   *reason = KSP_CONVERGED_ITERATING;
607: 
608:   KSPGetNormType(ksp,&normtype);
609:   if (normtype == KSP_NORM_NO) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Use KSPSkipConverged() with KSPNormType of KSP_NORM_NO");

611:   if (!ctx) SETERRQ(PETSC_ERR_ARG_NULL,"Convergence context must have been created with KSPDefaultConvergedCreate()");
612:   if (!n) {
613:     /* if user gives initial guess need to compute norm of b */
614:     if (!ksp->guess_zero && !ctx->initialrtol) {
615:       PetscReal      snorm;
616:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
617:         PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of RHS\n");
618:         VecNorm(ksp->vec_rhs,NORM_2,&snorm);        /*     <- b'*b */
619:       } else {
620:         Vec z;
621:         if (!ctx->work) {
622:           VecDuplicate(ksp->vec_rhs,&ctx->work);
623:         }
624:         z = ctx->work;
625:         KSP_PCApply(ksp,ksp->vec_rhs,z);
626:         if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
627:           PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n");
628:           VecNorm(z,NORM_2,&snorm);                 /*    dp <- b'*B'*B*b */
629:         } else if (ksp->normtype == KSP_NORM_NATURAL) {
630:           PetscScalar norm;
631:            PetscInfo(ksp,"user has provided nonzero initial guess, computing natural norm of RHS\n");
632:           VecDot(ksp->vec_rhs,z,&norm);
633:           snorm = sqrt(PetscAbsScalar(norm));                            /*    dp <- b'*B*b */
634:         }
635:       }
636:       /* handle special case of zero RHS and nonzero guess */
637:       if (!snorm) {
638:         PetscInfo(ksp,"Special case, user has provided nonzero initial guess and zero RHS\n");
639:         snorm = rnorm;
640:       }
641:       if (ctx->mininitialrtol) {
642:         ksp->rnorm0 = PetscMin(snorm,rnorm);
643:       } else {
644:         ksp->rnorm0 = snorm;
645:       }
646:     } else {
647:       ksp->rnorm0 = rnorm;
648:     }
649:     ksp->ttol   = PetscMax(ksp->rtol*ksp->rnorm0,ksp->abstol);
650:   }

652:   if (n <= ksp->chknorm) return(0);

654:   if (PetscIsInfOrNanScalar(rnorm)) {
655:     PetscInfo(ksp,"Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n");
656:     *reason = KSP_DIVERGED_NAN;
657:   } else if (rnorm <= ksp->ttol) {
658:     if (rnorm < ksp->abstol) {
659:       PetscInfo3(ksp,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G at iteration %D\n",rnorm,ksp->abstol,n);
660:       *reason = KSP_CONVERGED_ATOL;
661:     } else {
662:       if (ctx->initialrtol) {
663:         PetscInfo4(ksp,"Linear solver has converged. Residual norm %G is less than relative tolerance %G times initial residual norm %G at iteration %D\n",rnorm,ksp->rtol,ksp->rnorm0,n);
664:       } else {
665:         PetscInfo4(ksp,"Linear solver has converged. Residual norm %G is less than relative tolerance %G times initial right hand side norm %G at iteration %D\n",rnorm,ksp->rtol,ksp->rnorm0,n);
666:       }
667:       *reason = KSP_CONVERGED_RTOL;
668:     }
669:   } else if (rnorm >= ksp->divtol*ksp->rnorm0) {
670:     PetscInfo3(ksp,"Linear solver is diverging. Initial right hand size norm %G, current residual norm %G at iteration %D\n",ksp->rnorm0,rnorm,n);
671:     *reason = KSP_DIVERGED_DTOL;
672:   }
673:   return(0);
674: }

678: /*@C
679:    KSPDefaultConvergedDestroy - Frees the space used by the KSPDefaultConverged() function context

681:    Collective on KSP

683:    Input Parameters:
684: .  ctx - convergence context 

686:    Level: intermediate

688: .keywords: KSP, default, convergence, residual

690: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(),
691:           KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm(), KSPDefaultConvergedCreate()
692: @*/
693: PetscErrorCode  KSPDefaultConvergedDestroy(void *ctx)
694: {
695:   PetscErrorCode         ierr;
696:   KSPDefaultConvergedCtx *cctx = (KSPDefaultConvergedCtx*) ctx;

699:   if (cctx->work) {VecDestroy(cctx->work);}
700:   PetscFree(cctx);
701:   return(0);
702: }

706: /*
707:    KSPDefaultBuildSolution - Default code to create/move the solution.

709:    Input Parameters:
710: +  ksp - iterative context
711: -  v   - pointer to the user's vector  

713:    Output Parameter:
714: .  V - pointer to a vector containing the solution

716:    Level: advanced

718: .keywords:  KSP, build, solution, default

720: .seealso: KSPGetSolution(), KSPDefaultBuildResidual()
721: */
722: PetscErrorCode KSPDefaultBuildSolution(KSP ksp,Vec v,Vec *V)
723: {
726:   if (ksp->pc_side == PC_RIGHT) {
727:     if (ksp->pc) {
728:       if (v) {KSP_PCApply(ksp,ksp->vec_sol,v); *V = v;}
729:       else {SETERRQ(PETSC_ERR_SUP,"Not working with right preconditioner");}
730:     } else {
731:       if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
732:       else { *V = ksp->vec_sol;}
733:     }
734:   } else if (ksp->pc_side == PC_SYMMETRIC) {
735:     if (ksp->pc) {
736:       if (ksp->transpose_solve) SETERRQ(PETSC_ERR_SUP,"Not working with symmetric preconditioner and transpose solve");
737:       if (v) {PCApplySymmetricRight(ksp->pc,ksp->vec_sol,v); *V = v;}
738:       else {SETERRQ(PETSC_ERR_SUP,"Not working with symmetric preconditioner");}
739:     } else  {
740:       if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
741:       else { *V = ksp->vec_sol;}
742:     }
743:   } else {
744:     if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
745:     else { *V = ksp->vec_sol; }
746:   }
747:   return(0);
748: }

752: /*
753:    KSPDefaultBuildResidual - Default code to compute the residual.

755:    Input Parameters:
756: .  ksp - iterative context
757: .  t   - pointer to temporary vector
758: .  v   - pointer to user vector  

760:    Output Parameter:
761: .  V - pointer to a vector containing the residual

763:    Level: advanced

765: .keywords:  KSP, build, residual, default

767: .seealso: KSPDefaultBuildSolution()
768: */
769: PetscErrorCode KSPDefaultBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
770: {
772:   MatStructure   pflag;
773:   Mat            Amat,Pmat;

776:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
777:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
778:   KSPBuildSolution(ksp,t,PETSC_NULL);
779:   KSP_MatMult(ksp,Amat,t,v);
780:   VecAYPX(v,-1.0,ksp->vec_rhs);
781:   *V = v;
782:   return(0);
783: }

787: /*@C
788:   KSPGetVecs - Gets a number of work vectors.

790:   Input Parameters:
791: + ksp  - iterative context
792: . rightn  - number of right work vectors
793: - leftn   - number of left work vectors to allocate

795:   Output Parameter:
796: +  right - the array of vectors created
797: -  left - the array of left vectors

799:    Note: The right vector has as many elements as the matrix has columns. The left
800:      vector has as many elements as the matrix has rows.

802:    Level: advanced

804: .seealso:   MatGetVecs()

806: @*/
807: PetscErrorCode KSPGetVecs(KSP ksp,PetscInt rightn, Vec **right,PetscInt leftn,Vec **left)
808: {
810:   Vec            vecr,vecl;

813:   if (rightn) {
814:     if (!right) SETERRQ(PETSC_ERR_ARG_INCOMP,"You asked for right vectors but did not pass a pointer to hold them");
815:     if (ksp->vec_sol) vecr = ksp->vec_sol;
816:     else {
817:       Mat pmat;
818:       if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
819:       PCGetOperators(ksp->pc,PETSC_NULL,&pmat,PETSC_NULL);
820:       MatGetVecs(pmat,&vecr,PETSC_NULL);
821:     }
822:     VecDuplicateVecs(vecr,rightn,right);
823:     if (!ksp->vec_sol) {
824:       VecDestroy(vecr);
825:     }
826:   }
827:   if (leftn) {
828:     if (!left) SETERRQ(PETSC_ERR_ARG_INCOMP,"You asked for left vectors but did not pass a pointer to hold them");
829:     if (ksp->vec_rhs) vecl = ksp->vec_rhs;
830:     else {
831:       Mat pmat;
832:       if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
833:       PCGetOperators(ksp->pc,PETSC_NULL,&pmat,PETSC_NULL);
834:       MatGetVecs(pmat,PETSC_NULL,&vecl);
835:     }
836:     VecDuplicateVecs(vecl,leftn,left);
837:     if (!ksp->vec_rhs) {
838:       VecDestroy(vecl);
839:     }
840:   }
841:   return(0);
842: }

846: /*
847:   KSPDefaultGetWork - Gets a number of work vectors.

849:   Input Parameters:
850: . ksp  - iterative context
851: . nw   - number of work vectors to allocate

853:   Notes:
854:   Call this only if no work vectors have been allocated 
855:  */
856: PetscErrorCode KSPDefaultGetWork(KSP ksp,PetscInt nw)
857: {

861:   if (ksp->work) {KSPDefaultFreeWork(ksp);}
862:   ksp->nwork = nw;
863:   KSPGetVecs(ksp,nw,&ksp->work,0,PETSC_NULL);
864:   PetscLogObjectParents(ksp,nw,ksp->work);
865:   return(0);
866: }

870: /*
871:   KSPDefaultDestroy - Destroys a iterative context variable for methods with
872:   no separate context.  Preferred calling sequence KSPDestroy().

874:   Input Parameter: 
875: . ksp - the iterative context
876: */
877: PetscErrorCode KSPDefaultDestroy(KSP ksp)
878: {

883:   /* free work vectors */
884:   KSPDefaultFreeWork(ksp);
885:   /* free private data structure */
886:   PetscFree(ksp->data);
887:   return(0);
888: }

892: /*@
893:    KSPGetConvergedReason - Gets the reason the KSP iteration was stopped.

895:    Not Collective

897:    Input Parameter:
898: .  ksp - the KSP context

900:    Output Parameter:
901: .  reason - negative value indicates diverged, positive value converged, see KSPConvergedReason

903:    Possible values for reason:
904: +  KSP_CONVERGED_RTOL (residual 2-norm decreased by a factor of rtol, from 2-norm of right hand side)
905: .  KSP_CONVERGED_ATOL (residual 2-norm less than abstol)
906: .  KSP_CONVERGED_ITS (used by the preonly preconditioner that always uses ONE iteration, or when the KSPSkipConverged() convergence 
907:            test routine is set.
908: .  KSP_CONVERGED_CG_NEG_CURVE
909: .  KSP_CONVERGED_CG_CONSTRAINED
910: .  KSP_CONVERGED_STEP_LENGTH
911: .  KSP_DIVERGED_ITS  (required more than its to reach convergence)
912: .  KSP_DIVERGED_DTOL (residual norm increased by a factor of divtol)
913: .  KSP_DIVERGED_NAN (residual norm became Not-a-number likely do to 0/0)
914: .  KSP_DIVERGED_BREAKDOWN (generic breakdown in method)
915: -  KSP_DIVERGED_BREAKDOWN_BICG (Initial residual is orthogonal to preconditioned initial
916:                                 residual. Try a different preconditioner, or a different initial Level.)
917:  
918:    See also manual page for each reason.

920:    guess: beginner

922:    Notes: Can only be called after the call the KSPSolve() is complete.

924:    Level: intermediate
925:  
926: .keywords: KSP, nonlinear, set, convergence, test

928: .seealso: KSPSetConvergenceTest(), KSPDefaultConverged(), KSPSetTolerances(), KSPConvergedReason
929: @*/
930: PetscErrorCode  KSPGetConvergedReason(KSP ksp,KSPConvergedReason *reason)
931: {
935:   *reason = ksp->reason;
936:   return(0);
937: }