Actual source code: lsqr_converged.c

petsc-3.9.1 2018-04-29
Report Typos and Errors
  1:  #include <petsc/private/kspimpl.h>

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

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

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

 18:   *reason = KSP_CONVERGED_ITERATING;
 19:   if (iter == 0) return(0);

 21:   KSPCheckNorm(solksp,rnorm);

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

 29:   KSPGetSolution(solksp, &x_sol);
 30:   VecNorm(x_sol, NORM_2, &xnorm);

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

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

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

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