Actual source code: iterativ.c

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

 12: /*@
 13:    KSPGetResidualNorm - Gets the last (approximate preconditioned)
 14:    residual norm that has been computed.
 15:  
 16:    Not Collective

 18:    Input Parameters:
 19: .  ksp - the iterative context

 21:    Output Parameters:
 22: .  rnorm - residual norm

 24:    Level: intermediate

 26: .keywords: KSP, get, residual norm

 28: .seealso: KSPBuildResidual()
 29: @*/
 30: PetscErrorCode  KSPGetResidualNorm(KSP ksp,PetscReal *rnorm)
 31: {
 35:   *rnorm = ksp->rnorm;
 36:   return(0);
 37: }

 41: /*@
 42:    KSPGetIterationNumber - Gets the current iteration number; if the 
 43:          KSPSolve() is complete, returns the number of iterations
 44:          used.
 45:  
 46:    Not Collective

 48:    Input Parameters:
 49: .  ksp - the iterative context

 51:    Output Parameters:
 52: .  its - number of iterations

 54:    Level: intermediate

 56:    Notes:
 57:       During the ith iteration this returns i-1
 58: .keywords: KSP, get, residual norm

 60: .seealso: KSPBuildResidual(), KSPGetResidualNorm()
 61: @*/
 62: PetscErrorCode  KSPGetIterationNumber(KSP ksp,PetscInt *its)
 63: {
 67:   *its = ksp->its;
 68:   return(0);
 69: }

 73: /*@C
 74:     KSPMonitorSingularValue - Prints the two norm of the true residual and
 75:     estimation of the extreme singular values of the preconditioned problem
 76:     at each iteration.
 77:  
 78:     Logically Collective on KSP

 80:     Input Parameters:
 81: +   ksp - the iterative context
 82: .   n  - the iteration
 83: -   rnorm - the two norm of the residual

 85:     Options Database Key:
 86: .   -ksp_monitor_singular_value - Activates KSPMonitorSingularValue()

 88:     Notes:
 89:     The CG solver uses the Lanczos technique for eigenvalue computation, 
 90:     while GMRES uses the Arnoldi technique; other iterative methods do
 91:     not currently compute singular values.

 93:     Level: intermediate

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

 97: .seealso: KSPComputeExtremeSingularValues()
 98: @*/
 99: PetscErrorCode  KSPMonitorSingularValue(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
100: {
101:   PetscReal       emin,emax,c;
102:   PetscErrorCode  ierr;
103:   PetscViewer     viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm);

107:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
108:   if (!ksp->calc_sings) {
109:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,(double)rnorm);
110:   } else {
111:     KSPComputeExtremeSingularValues(ksp,&emax,&emin);
112:     c = emax/emin;
113:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e %% max %14.12e min %14.12e max/min %14.12e\n",n,(double)rnorm,(double)emax,(double)emin,(double)c);
114:   }
115:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
116:   return(0);
117: }

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

125:    Collective on KSP

127:    Input Parameters:
128: +  ksp - the KSP context
129: .  its - iteration number
130: .  fgnorm - 2-norm of residual (or gradient)
131: -  dummy - either a viewer or PETSC_NULL

133:    Level: intermediate

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

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

141: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView()
142: @*/
143: PetscErrorCode  KSPMonitorSolution(KSP ksp,PetscInt its,PetscReal fgnorm,void *dummy)
144: {
146:   Vec            x;
147:   PetscViewer    viewer = (PetscViewer) dummy;

150:   KSPBuildSolution(ksp,PETSC_NULL,&x);
151:   if (!viewer) {
152:     MPI_Comm comm;
153:     PetscObjectGetComm((PetscObject)ksp,&comm);
154:     viewer = PETSC_VIEWER_DRAW_(comm);
155:   }
156:   VecView(x,viewer);

158:   return(0);
159: }

163: /*@C
164:    KSPMonitorDefault - Print the residual norm at each iteration of an
165:    iterative solver.

167:    Collective on KSP

169:    Input Parameters:
170: +  ksp   - iterative context
171: .  n     - iteration number
172: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
173: -  dummy - unused monitor context 

175:    Level: intermediate

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

179: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGCreate()
180: @*/
181: PetscErrorCode  KSPMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
182: {
184:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm);

187:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
188:   if (n == 0 && ((PetscObject)ksp)->prefix) {
189:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
190:   }
191:   PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,(double)rnorm);
192:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
193:   return(0);
194: }

198: /*@C
199:    KSPMonitorTrueResidualNorm - Prints the true residual norm as well as the preconditioned
200:    residual norm at each iteration of an iterative solver.

202:    Collective on KSP

204:    Input Parameters:
205: +  ksp   - iterative context
206: .  n     - iteration number
207: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
208: -  dummy - unused monitor context 

210:    Options Database Key:
211: .  -ksp_monitor_true_residual - Activates KSPMonitorTrueResidualNorm()

213:    Notes:
214:    When using right preconditioning, these values are equivalent.

216:    Level: intermediate

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

220: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGCreate()
221: @*/
222: PetscErrorCode  KSPMonitorTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
223: {
224:   PetscErrorCode  ierr;
225:   Vec             resid,work;
226:   PetscReal       scnorm,bnorm;
227:   PC              pc;
228:   Mat             A,B;
229:   PetscViewer     viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm);
230:   char            normtype[256];

233:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
234:   if (n == 0 && ((PetscObject)ksp)->prefix) {
235:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
236:   }
237:   VecDuplicate(ksp->vec_rhs,&work);
238:   KSPBuildResidual(ksp,0,work,&resid);

240:   /*
241:      Unscale the residual but only if both matrices are the same matrix, since only then would 
242:     they be scaled.
243:   */
244:   VecCopy(resid,work);
245:   KSPGetPC(ksp,&pc);
246:   PCGetOperators(pc,&A,&B,PETSC_NULL);
247:   if (A == B) {
248:     MatUnScaleSystem(A,work,PETSC_NULL);
249:   }
250:   VecNorm(work,NORM_2,&scnorm);
251:   VecDestroy(&work);
252:   VecNorm(ksp->vec_rhs,NORM_2,&bnorm);
253:   PetscStrncpy(normtype,KSPNormTypes[ksp->normtype],sizeof normtype);
254:   PetscStrtolower(normtype);
255:   PetscViewerASCIIPrintf(viewer,"%3D KSP %s resid norm %14.12e true resid norm %14.12e ||r(i)||/||b|| %14.12e\n",n,normtype,(double)rnorm,(double)scnorm,(double)(scnorm/bnorm));
256:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
257:   return(0);
258: }

262: PetscErrorCode  KSPMonitorRange_Private(KSP ksp,PetscInt it,PetscReal *per)
263: {
264:   PetscErrorCode          ierr;
265:   Vec                     resid,work;
266:   PetscReal               rmax,pwork;
267:   PC                      pc;
268:   Mat                     A,B;
269:   PetscInt                i,n,N;
270:   PetscScalar             *r;

273:   VecDuplicate(ksp->vec_rhs,&work);
274:   KSPBuildResidual(ksp,0,work,&resid);

276:   /*
277:      Unscale the residual if the matrix is, but only if both matrices are the same matrix, since only then would 
278:     they be scaled.
279:   */
280:   VecCopy(resid,work);
281:   KSPGetPC(ksp,&pc);
282:   PCGetOperators(pc,&A,&B,PETSC_NULL);
283:   if (A == B) {
284:     MatUnScaleSystem(A,work,PETSC_NULL);
285:   }
286:   VecNorm(work,NORM_INFINITY,&rmax);
287:   VecGetLocalSize(work,&n);
288:   VecGetSize(work,&N);
289:   VecGetArray(work,&r);
290:   pwork = 0.0;
291:   for (i=0; i<n; i++) {
292:     pwork += (PetscAbsScalar(r[i]) > .20*rmax);
293:   }
294:   MPI_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,((PetscObject)ksp)->comm);
295:   VecRestoreArray(work,&r);
296:   VecDestroy(&work);
297:   *per  = *per/N;
298:   return(0);
299: }

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

306:    Collective on KSP

308:    Input Parameters:
309: +  ksp   - iterative context
310: .  it    - iteration number
311: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).  
312: -  dummy - unused monitor context 

314:    Options Database Key:
315: .  -ksp_monitor_range - Activates KSPMonitorRange()

317:    Level: intermediate

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

321: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGCreate()
322: @*/
323: PetscErrorCode  KSPMonitorRange(KSP ksp,PetscInt it,PetscReal rnorm,void *dummy)
324: {
325:   PetscErrorCode   ierr;
326:   PetscReal        perc,rel;
327:   PetscViewer      viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm);
328:   /* should be in a MonitorRangeContext */
329:   static PetscReal prev;

332:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
333:   if (!it) prev = rnorm;
334:   if (it == 0 && ((PetscObject)ksp)->prefix) {
335:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
336:   }
337:   KSPMonitorRange_Private(ksp,it,&perc);

339:   rel  = (prev - rnorm)/prev;
340:   prev = rnorm;
341:   PetscViewerASCIIPrintf(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,(double)rnorm,(double)100.0*perc,(double)rel,(double)rel/perc);
342:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
343:   return(0);
344: }

348: /*
349:   Default (short) KSP Monitor, same as KSPMonitorDefault() except
350:   it prints fewer digits of the residual as the residual gets smaller.
351:   This is because the later digits are meaningless and are often 
352:   different on different machines; by using this routine different 
353:   machines will usually generate the same output.
354: */
355: PetscErrorCode  KSPMonitorDefaultShort(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy)
356: {
358:   PetscViewer    viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(((PetscObject)ksp)->comm);

361:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
362:   if (its == 0 && ((PetscObject)ksp)->prefix) {
363:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
364:   }

366:   if (fnorm > 1.e-9) {
367:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %G \n",its,fnorm);
368:   } else if (fnorm > 1.e-11){
369:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %5.3e \n",its,(double)fnorm);
370:   } else {
371:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm < 1.e-11\n",its);
372:   }
373:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
374:   return(0);
375: }

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

383:    Collective on KSP

385:    Input Parameters:
386: +  ksp   - iterative context
387: .  n     - iteration number
388: .  rnorm - 2-norm residual value (may be estimated)
389: -  dummy - unused convergence context 

391:    Returns:
392: .  reason - KSP_CONVERGED_ITERATING, KSP_CONVERGED_ITS

394:    Notes: 
395:    This should be used as the convergence test with the option
396:    KSPSetNormType(ksp,KSP_NORM_NONE), since norms of the residual are
397:    not computed. Convergence is then declared after the maximum number
398:    of iterations have been reached. Useful when one is using CG or
399:    BiCGStab as a smoother.
400:                     
401:    Level: advanced

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

405: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSetNormType()
406: @*/
407: PetscErrorCode  KSPSkipConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
408: {
412:   *reason = KSP_CONVERGED_ITERATING;
413:   if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
414:   return(0);
415: }


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

423:    Collective on KSP

425:    Output Parameter:
426: .  ctx - convergence context 

428:    Level: intermediate

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

432: .seealso: KSPDefaultConverged(), KSPDefaultConvergedDestroy(), KSPSetConvergenceTest(), KSPSetTolerances(),
433:           KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm()
434: @*/
435: PetscErrorCode  KSPDefaultConvergedCreate(void **ctx)
436: {
437:   PetscErrorCode         ierr;
438:   KSPDefaultConvergedCtx *cctx;

441:   PetscNew(KSPDefaultConvergedCtx,&cctx);
442:   *ctx = cctx;
443:   return(0);
444: }

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

453:    Collective on KSP

455:    Input Parameters:
456: .  ksp   - iterative context

458:    Options Database:
459: .   -ksp_converged_use_initial_residual_norm

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

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

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

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

471: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUMIRNorm()
472: @*/
473: PetscErrorCode  KSPDefaultConvergedSetUIRNorm(KSP ksp)
474: {
475:   KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) ksp->cnvP;

479:   if (ksp->converged != KSPDefaultConverged) return(0);
480:   if (ctx->mininitialrtol) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPDefaultConvergedSetUIRNorm() and KSPDefaultConvergedSetUMIRNorm() together");
481:   ctx->initialrtol = PETSC_TRUE;
482:   return(0);
483: }

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

492:    Collective on KSP

494:    Input Parameters:
495: .  ksp   - iterative context

497:    Options Database:
498: .   -ksp_converged_use_min_initial_residual_norm

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

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

505:    Level: intermediate

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

509: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm()
510: @*/
511: PetscErrorCode  KSPDefaultConvergedSetUMIRNorm(KSP ksp)
512: {
513:   KSPDefaultConvergedCtx *ctx = (KSPDefaultConvergedCtx*) ksp->cnvP;

517:   if (ksp->converged != KSPDefaultConverged) return(0);
518:   if (ctx->initialrtol) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPDefaultConvergedSetUIRNorm() and KSPDefaultConvergedSetUMIRNorm() together");
519:   ctx->mininitialrtol = PETSC_TRUE;
520:    return(0);
521: }

525: /*@C
526:    KSPDefaultConverged - Determines convergence of
527:    the iterative solvers (default code).

529:    Collective on KSP

531:    Input Parameters:
532: +  ksp   - iterative context
533: .  n     - iteration number
534: .  rnorm - 2-norm residual value (may be estimated)
535: -  ctx - convergence context which must be created by KSPDefaultConvergedCreate()

537:    reason is set to:
538: +   positive - if the iteration has converged;
539: .   negative - if residual norm exceeds divergence threshold;
540: -   0 - otherwise.

542:    Notes:
543:    KSPDefaultConverged() reaches convergence when
544: $      rnorm < MAX (rtol * rnorm_0, abstol);
545:    Divergence is detected if
546: $      rnorm > dtol * rnorm_0,

548:    where 
549: +     rtol = relative tolerance,
550: .     abstol = absolute tolerance.
551: .     dtol = divergence tolerance,
552: -     rnorm_0 is the two norm of the right hand side. When initial guess is non-zero you 
553:           can call KSPDefaultConvergedSetUIRNorm() to use the norm of (b - A*(initial guess))
554:           as the starting point for relative norm convergence testing.

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

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

561:    Level: intermediate

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

565: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(), KSPConvergedReason, KSPGetConvergedReason(),
566:           KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm(), KSPDefaultConvergedCreate(), KSPDefaultConvergedDestroy()
567: @*/
568: PetscErrorCode  KSPDefaultConverged(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
569: {
570:   PetscErrorCode         ierr;
571:   KSPDefaultConvergedCtx *cctx = (KSPDefaultConvergedCtx*) ctx;
572:   KSPNormType            normtype;

577:   *reason = KSP_CONVERGED_ITERATING;
578: 
579:   KSPGetNormType(ksp,&normtype);
580:   if (normtype == KSP_NORM_NONE) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_WRONGSTATE,"Use KSPSkipConverged() with KSPNormType of KSP_NORM_NONE");

582:   if (!cctx) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_NULL,"Convergence context must have been created with KSPDefaultConvergedCreate()");
583:   if (!n) {
584:     /* if user gives initial guess need to compute norm of b */
585:     if (!ksp->guess_zero && !cctx->initialrtol) {
586:       PetscReal      snorm;
587:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
588:         PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of RHS\n");
589:         VecNorm(ksp->vec_rhs,NORM_2,&snorm);        /*     <- b'*b */
590:       } else {
591:         Vec z;
592:         /* Should avoid allocating the z vector each time but cannot stash it in cctx because if KSPReset() is called the vector size might change */
593:         VecDuplicate(ksp->vec_rhs,&z);
594:         KSP_PCApply(ksp,ksp->vec_rhs,z);
595:         if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
596:           PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n");
597:           VecNorm(z,NORM_2,&snorm);                 /*    dp <- b'*B'*B*b */
598:         } else if (ksp->normtype == KSP_NORM_NATURAL) {
599:           PetscScalar norm;
600:            PetscInfo(ksp,"user has provided nonzero initial guess, computing natural norm of RHS\n");
601:           VecDot(ksp->vec_rhs,z,&norm);
602:           snorm = PetscSqrtReal(PetscAbsScalar(norm));                            /*    dp <- b'*B*b */
603:         }
604:         VecDestroy(&z);
605:       }
606:       /* handle special case of zero RHS and nonzero guess */
607:       if (!snorm) {
608:         PetscInfo(ksp,"Special case, user has provided nonzero initial guess and zero RHS\n");
609:         snorm = rnorm;
610:       }
611:       if (cctx->mininitialrtol) {
612:         ksp->rnorm0 = PetscMin(snorm,rnorm);
613:       } else {
614:         ksp->rnorm0 = snorm;
615:       }
616:     } else {
617:       ksp->rnorm0 = rnorm;
618:     }
619:     ksp->ttol   = PetscMax(ksp->rtol*ksp->rnorm0,ksp->abstol);
620:   }

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

624:   if (PetscIsInfOrNanScalar(rnorm)) {
625:     PetscInfo(ksp,"Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n");
626:     *reason = KSP_DIVERGED_NAN;
627:   } else if (rnorm <= ksp->ttol) {
628:     if (rnorm < ksp->abstol) {
629:       PetscInfo3(ksp,"Linear solver has converged. Residual norm %14.12e is less than absolute tolerance %14.12e at iteration %D\n",(double)rnorm,(double)ksp->abstol,n);
630:       *reason = KSP_CONVERGED_ATOL;
631:     } else {
632:       if (cctx->initialrtol) {
633:         PetscInfo4(ksp,"Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial residual norm %14.12e at iteration %D\n",(double)rnorm,(double)ksp->rtol,(double)ksp->rnorm0,n);
634:       } else {
635:         PetscInfo4(ksp,"Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial right hand side norm %14.12e at iteration %D\n",(double)rnorm,(double)ksp->rtol,(double)ksp->rnorm0,n);
636:       }
637:       *reason = KSP_CONVERGED_RTOL;
638:     }
639:   } else if (rnorm >= ksp->divtol*ksp->rnorm0) {
640:     PetscInfo3(ksp,"Linear solver is diverging. Initial right hand size norm %14.12e, current residual norm %14.12e at iteration %D\n",(double)ksp->rnorm0,(double)rnorm,n);
641:     *reason = KSP_DIVERGED_DTOL;
642:   }
643:   return(0);
644: }

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

651:    Collective on KSP

653:    Input Parameters:
654: .  ctx - convergence context 

656:    Level: intermediate

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

660: .seealso: KSPDefaultConverged(), KSPDefaultConvergedCreate(), KSPSetConvergenceTest(), KSPSetTolerances(), KSPSkipConverged(),
661:           KSPConvergedReason, KSPGetConvergedReason(), KSPDefaultConvergedSetUIRNorm(), KSPDefaultConvergedSetUMIRNorm()
662: @*/
663: PetscErrorCode  KSPDefaultConvergedDestroy(void *ctx)
664: {
665:   PetscErrorCode         ierr;
666:   KSPDefaultConvergedCtx *cctx = (KSPDefaultConvergedCtx*) ctx;

669:   VecDestroy(&cctx->work);
670:   PetscFree(ctx);
671:   return(0);
672: }

676: /*
677:    KSPDefaultBuildSolution - Default code to create/move the solution.

679:    Input Parameters:
680: +  ksp - iterative context
681: -  v   - pointer to the user's vector  

683:    Output Parameter:
684: .  V - pointer to a vector containing the solution

686:    Level: advanced

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

690: .seealso: KSPGetSolution(), KSPDefaultBuildResidual()
691: */
692: PetscErrorCode KSPDefaultBuildSolution(KSP ksp,Vec v,Vec *V)
693: {
696:   if (ksp->pc_side == PC_RIGHT) {
697:     if (ksp->pc) {
698:       if (v) {KSP_PCApply(ksp,ksp->vec_sol,v); *V = v;}
699:       else SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Not working with right preconditioner");
700:     } else {
701:       if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
702:       else { *V = ksp->vec_sol;}
703:     }
704:   } else if (ksp->pc_side == PC_SYMMETRIC) {
705:     if (ksp->pc) {
706:       if (ksp->transpose_solve) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Not working with symmetric preconditioner and transpose solve");
707:       if (v) {PCApplySymmetricRight(ksp->pc,ksp->vec_sol,v); *V = v;}
708:       else SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Not working with symmetric preconditioner");
709:     } else  {
710:       if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
711:       else { *V = ksp->vec_sol;}
712:     }
713:   } else {
714:     if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
715:     else { *V = ksp->vec_sol; }
716:   }
717:   return(0);
718: }

722: /*
723:    KSPDefaultBuildResidual - Default code to compute the residual.

725:    Input Parameters:
726: .  ksp - iterative context
727: .  t   - pointer to temporary vector
728: .  v   - pointer to user vector  

730:    Output Parameter:
731: .  V - pointer to a vector containing the residual

733:    Level: advanced

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

737: .seealso: KSPDefaultBuildSolution()
738: */
739: PetscErrorCode KSPDefaultBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
740: {
742:   MatStructure   pflag;
743:   Mat            Amat,Pmat;

746:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
747:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
748:   KSPBuildSolution(ksp,t,PETSC_NULL);
749:   KSP_MatMult(ksp,Amat,t,v);
750:   VecAYPX(v,-1.0,ksp->vec_rhs);
751:   *V = v;
752:   return(0);
753: }

757: /*@C
758:   KSPGetVecs - Gets a number of work vectors.

760:   Input Parameters:
761: + ksp  - iterative context
762: . rightn  - number of right work vectors
763: - leftn   - number of left work vectors to allocate

765:   Output Parameter:
766: +  right - the array of vectors created
767: -  left - the array of left vectors

769:    Note: The right vector has as many elements as the matrix has columns. The left
770:      vector has as many elements as the matrix has rows.

772:    Level: advanced

774: .seealso:   MatGetVecs()

776: @*/
777: PetscErrorCode KSPGetVecs(KSP ksp,PetscInt rightn, Vec **right,PetscInt leftn,Vec **left)
778: {
780:   Vec            vecr,vecl;

783:   if (rightn) {
784:     if (!right) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_INCOMP,"You asked for right vectors but did not pass a pointer to hold them");
785:     if (ksp->vec_sol) vecr = ksp->vec_sol;
786:     else {
787:       if (ksp->dm) {
788:         DMGetGlobalVector(ksp->dm,&vecr);
789:       } else {
790:         Mat pmat;
791:         if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
792:         PCGetOperators(ksp->pc,PETSC_NULL,&pmat,PETSC_NULL);
793:         MatGetVecs(pmat,&vecr,PETSC_NULL);
794:       }
795:     }
796:     VecDuplicateVecs(vecr,rightn,right);
797:     if (!ksp->vec_sol) {
798:       if (ksp->dm) {
799:         DMRestoreGlobalVector(ksp->dm,&vecr);
800:       } else {
801:         VecDestroy(&vecr);
802:       }
803:     }
804:   }
805:   if (leftn) {
806:     if (!left) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_ARG_INCOMP,"You asked for left vectors but did not pass a pointer to hold them");
807:     if (ksp->vec_rhs) vecl = ksp->vec_rhs;
808:     else {
809:       if (ksp->dm) {
810:         DMGetGlobalVector(ksp->dm,&vecl);
811:       } else {
812:         Mat pmat;
813:         if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
814:         PCGetOperators(ksp->pc,PETSC_NULL,&pmat,PETSC_NULL);
815:         MatGetVecs(pmat,PETSC_NULL,&vecl);
816:       }
817:     }
818:     VecDuplicateVecs(vecl,leftn,left);
819:     if (!ksp->vec_rhs) {
820:       if (ksp->dm) {
821:         DMRestoreGlobalVector(ksp->dm,&vecl);
822:       } else {
823:         VecDestroy(&vecl);
824:       }
825:     }
826:   }
827:   return(0);
828: }

832: /*
833:   KSPDefaultGetWork - Gets a number of work vectors.

835:   Input Parameters:
836: . ksp  - iterative context
837: . nw   - number of work vectors to allocate

839:   Notes:
840:   Call this only if no work vectors have been allocated 
841:  */
842: PetscErrorCode KSPDefaultGetWork(KSP ksp,PetscInt nw)
843: {

847:   VecDestroyVecs(ksp->nwork,&ksp->work);
848:   ksp->nwork = nw;
849:   KSPGetVecs(ksp,nw,&ksp->work,0,PETSC_NULL);
850:   PetscLogObjectParents(ksp,nw,ksp->work);
851:   return(0);
852: }

856: /*
857:   KSPDefaultDestroy - Destroys a iterative context variable for methods with
858:   no separate context.  Preferred calling sequence KSPDestroy().

860:   Input Parameter:
861: . ksp - the iterative context
862: */
863: PetscErrorCode KSPDefaultDestroy(KSP ksp)
864: {

869:   PetscFree(ksp->data);
870:   return(0);
871: }

875: /*@
876:    KSPGetConvergedReason - Gets the reason the KSP iteration was stopped.

878:    Not Collective

880:    Input Parameter:
881: .  ksp - the KSP context

883:    Output Parameter:
884: .  reason - negative value indicates diverged, positive value converged, see KSPConvergedReason

886:    Possible values for reason:
887: +  KSP_CONVERGED_RTOL (residual 2-norm decreased by a factor of rtol, from 2-norm of right hand side)
888: .  KSP_CONVERGED_ATOL (residual 2-norm less than abstol)
889: .  KSP_CONVERGED_ITS (used by the preonly preconditioner that always uses ONE iteration, or when the KSPSkipConverged() convergence 
890:            test routine is set.
891: .  KSP_CONVERGED_CG_NEG_CURVE
892: .  KSP_CONVERGED_CG_CONSTRAINED
893: .  KSP_CONVERGED_STEP_LENGTH
894: .  KSP_DIVERGED_ITS  (required more than its to reach convergence)
895: .  KSP_DIVERGED_DTOL (residual norm increased by a factor of divtol)
896: .  KSP_DIVERGED_NAN (residual norm became Not-a-number likely due to 0/0)
897: .  KSP_DIVERGED_BREAKDOWN (generic breakdown in method)
898: -  KSP_DIVERGED_BREAKDOWN_BICG (Initial residual is orthogonal to preconditioned initial
899:                                 residual. Try a different preconditioner, or a different initial Level.)
900:  
901:    See also manual page for each reason.

903:    guess: beginner

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

907:    Level: intermediate
908:  
909: .keywords: KSP, nonlinear, set, convergence, test

911: .seealso: KSPSetConvergenceTest(), KSPDefaultConverged(), KSPSetTolerances(), KSPConvergedReason
912: @*/
913: PetscErrorCode  KSPGetConvergedReason(KSP ksp,KSPConvergedReason *reason)
914: {
918:   *reason = ksp->reason;
919:   return(0);
920: }

924: /*@
925:    KSPSetDM - Sets the DM that may be used by some preconditioners

927:    Logically Collective on KSP

929:    Input Parameters:
930: +  ksp - the preconditioner context
931: -  dm - the dm

933:    Level: intermediate


936: .seealso: KSPGetDM(), KSPSetDM(), KSPGetDM()
937: @*/
938: PetscErrorCode  KSPSetDM(KSP ksp,DM dm)
939: {
941:   PC             pc;

945:   if (dm) {PetscObjectReference((PetscObject)dm);}
946:   DMDestroy(&ksp->dm);
947:   ksp->dm = dm;
948:   KSPGetPC(ksp,&pc);
949:   PCSetDM(pc,dm);
950:   ksp->dmActive = PETSC_TRUE;
951:   return(0);
952: }

956: /*@
957:    KSPSetDMActive - Indicates the DM should be used to generate the linear system matrix and right hand side

959:    Logically Collective on KSP

961:    Input Parameters:
962: +  ksp - the preconditioner context
963: -  flg - use the DM

965:    Level: intermediate

967:    Notes:
968:    By default KSPSetDM() sets the DM as active, call KSPSetDMActive(dm,PETSC_FALSE); after KSPSetDM(dm) to not have the KSP object use the DM to generate the matrices

970: .seealso: KSPGetDM(), KSPSetDM(), KSPGetDM()
971: @*/
972: PetscErrorCode  KSPSetDMActive(KSP ksp,PetscBool  flg)
973: {
977:   ksp->dmActive = flg;
978:   return(0);
979: }

983: /*@
984:    KSPGetDM - Gets the DM that may be used by some preconditioners

986:    Not Collective

988:    Input Parameter:
989: . ksp - the preconditioner context

991:    Output Parameter:
992: .  dm - the dm

994:    Level: intermediate


997: .seealso: KSPSetDM(), KSPSetDM(), KSPGetDM()
998: @*/
999: PetscErrorCode  KSPGetDM(KSP ksp,DM *dm)
1000: {
1003:   *dm = ksp->dm;
1004:   return(0);
1005: }

1009: /*@
1010:    KSPSetApplicationContext - Sets the optional user-defined context for the linear solver.

1012:    Logically Collective on KSP

1014:    Input Parameters:
1015: +  ksp - the KSP context
1016: -  usrP - optional user context

1018:    Level: intermediate

1020: .keywords: KSP, set, application, context

1022: .seealso: KSPGetApplicationContext()
1023: @*/
1024: PetscErrorCode  KSPSetApplicationContext(KSP ksp,void *usrP)
1025: {
1027:   PC             pc;

1031:   ksp->user = usrP;
1032:   KSPGetPC(ksp,&pc);
1033:   PCSetApplicationContext(pc,usrP);
1034:   return(0);
1035: }

1039: /*@
1040:    KSPGetApplicationContext - Gets the user-defined context for the linear solver.

1042:    Not Collective

1044:    Input Parameter:
1045: .  ksp - KSP context

1047:    Output Parameter:
1048: .  usrP - user context

1050:    Level: intermediate

1052: .keywords: KSP, get, application, context

1054: .seealso: KSPSetApplicationContext()
1055: @*/
1056: PetscErrorCode  KSPGetApplicationContext(KSP ksp,void *usrP)
1057: {
1060:   *(void**)usrP = ksp->user;
1061:   return(0);
1062: }