Actual source code: lsqr_converged.c

petsc-3.7.3 2016-07-24
Report Typos and Errors
  1: #include <petsc/private/kspimpl.h>
  2: #include <../src/ksp/ksp/impls/lsqr/lsqr.h>
  3: extern PetscErrorCode  KSPLSQRGetArnorm(KSP,PetscReal*,PetscReal*,PetscReal*);

  5: PetscErrorCode KSPConvergedLSQR(KSP solksp,PetscInt iter,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
  6: {
  7:   PetscReal atol;          /* Absolute convergence criterium */
  8:   PetscReal dtol;          /* Divergence criterium */
  9:   PetscReal rtol;          /* Relative convergence criterium */
 10:   PetscReal stop1;         /* Stop if: |r| < rtol*|b| + atol*|A|*|x| */
 11:   PetscReal stop2;         /* Stop if: |A^t.r|/(|A|*|r|) < atol */
 12:   Vec       x_sol;         /* Current solution vector */

 14:   PetscReal arnorm, anorm, bnorm, xnorm;        /* Norms of A*residual; matrix A; rhs; solution */

 16:   PetscInt       mxiter;   /* Maximum # of iterations */

 20:   *reason = KSP_CONVERGED_ITERATING;
 21:   if (iter == 0) return(0);

 23:   KSPCheckNorm(solksp,rnorm);

 25:   KSPGetTolerances(solksp, &rtol, &atol, &dtol, &mxiter);
 26:   if (iter > mxiter) {
 27:     *reason = KSP_DIVERGED_ITS;
 28:     return(0);
 29:   }

 31:   KSPGetSolution(solksp, &x_sol);
 32:   VecNorm(x_sol, NORM_2, &xnorm);

 34:   KSPLSQRGetArnorm(solksp, &arnorm, &bnorm, &anorm);
 35:   if (bnorm > 0.0) {
 36:     stop1 = rnorm / bnorm;
 37:     rtol  = rtol + atol * anorm*xnorm/bnorm;
 38:   } else {
 39:     stop1 = 0.0;
 40:     rtol  = 0.0;
 41:   }
 42:   stop2 = 0.0;
 43:   if (rnorm > 0.0) stop2 = arnorm / (anorm * rnorm);

 45:   /* Test for tolerances set by the user */
 46:   if (stop1 <= rtol) *reason = KSP_CONVERGED_RTOL;
 47:   if (stop2 <= atol) *reason = KSP_CONVERGED_ATOL;

 49:   /* Test for machine precision */
 50:   if (bnorm > 0) stop1 = stop1 / (1.0 + anorm*xnorm/bnorm);
 51:   else stop1 = 0.0;

 53:   stop1 = 1.0 + stop1;
 54:   stop2 = 1.0 + stop2;
 55:   if (stop1 <= 1.0) *reason = KSP_CONVERGED_RTOL;
 56:   if (stop2 <= 1.0) *reason = KSP_CONVERGED_ATOL;
 57:   return(0);
 58: }