Actual source code: bcgs.c

  1: /*$Id: bcgs.c,v 1.78 2001/08/07 03:03:49 balay Exp $*/

  3: /*                       
  4:     This code implements the BiCGStab (Stabilized version of BiConjugate
  5:     Gradient Squared) method.  Reference: van der Vorst, SIAM J. Sci. Stat. Comput., 1992.

  7:     Note that for the complex numbers version, the VecDot() arguments
  8:     within the code MUST remain in the order given for correct computation
  9:     of inner products.
 10: */
 11:  #include src/sles/ksp/kspimpl.h

 13: static int KSPSetUp_BCGS(KSP ksp)
 14: {

 18:   if (ksp->pc_side == PC_SYMMETRIC) {
 19:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPBCGS");
 20:   }
 21:   KSPDefaultGetWork(ksp,6);
 22:   return(0);
 23: }

 25: static int  KSPSolve_BCGS(KSP ksp,int *its)
 26: {
 27:   int         i,maxit,ierr;
 28:   PetscScalar rho,rhoold,alpha,beta,omega,omegaold,d1,d2,zero = 0.0,tmp;
 29:   Vec         X,B,V,P,R,RP,T,S;
 30:   PetscReal   dp = 0.0;


 34:   maxit   = ksp->max_it;
 35:   X       = ksp->vec_sol;
 36:   B       = ksp->vec_rhs;
 37:   R       = ksp->work[0];
 38:   RP      = ksp->work[1];
 39:   V       = ksp->work[2];
 40:   T       = ksp->work[3];
 41:   S       = ksp->work[4];
 42:   P       = ksp->work[5];

 44:   /* Compute initial preconditioned residual */
 45:   KSPInitialResidual(ksp,X,V,T,R,B);

 47:   /* Test for nothing to do */
 48:   if (ksp->normtype != KSP_NO_NORM) {
 49:     VecNorm(R,NORM_2,&dp);
 50:   }
 51:   PetscObjectTakeAccess(ksp);
 52:   ksp->its   = 0;
 53:   ksp->rnorm = dp;
 54:   PetscObjectGrantAccess(ksp);
 55:   KSPLogResidualHistory(ksp,dp);
 56:   KSPMonitor(ksp,0,dp);
 57:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 58:   if (ksp->reason) {*its = 0; return(0);}

 60:   /* Make the initial Rp == R */
 61:   VecCopy(R,RP);

 63:   rhoold   = 1.0;
 64:   alpha    = 1.0;
 65:   omegaold = 1.0;
 66:   VecSet(&zero,P);
 67:   VecSet(&zero,V);

 69:   for (i=0; i<maxit; i++) {

 71:     VecDot(R,RP,&rho);       /*   rho <- (r,rp)      */
 72:     if (rho == 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Breakdown, rho = r . rp = 0");
 73:     beta = (rho/rhoold) * (alpha/omegaold);
 74:     tmp = -omegaold; VecAXPY(&tmp,V,P);            /*   p <- p - w v       */
 75:     VecAYPX(&beta,R,P);      /*   p <- r + p beta    */
 76:     KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T);  /*   v <- K p           */
 77:     VecDot(V,RP,&d1);
 78:     alpha = rho / d1; tmp = -alpha;                /*   a <- rho / (v,rp)  */
 79:     VecWAXPY(&tmp,V,R,S);    /*   s <- r - a v       */
 80:     KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,S,T,R);/*   t <- K s    */
 81:     VecDot(S,T,&d1);
 82:     VecDot(T,T,&d2);
 83:     if (d2 == 0.0) {
 84:       /* t is 0.  if s is 0, then alpha v == r, and hence alpha p
 85:          may be our solution.  Give it a try? */
 86:       VecDot(S,S,&d1);
 87:       if (d1 != 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Breakdown, da = s . s != 0");
 88:       VecAXPY(&alpha,P,X);   /*   x <- x + a p       */
 89:       PetscObjectTakeAccess(ksp);
 90:       ksp->its++;
 91:       ksp->rnorm  = 0.0;
 92:       ksp->reason = KSP_CONVERGED_RTOL;
 93:       PetscObjectGrantAccess(ksp);
 94:       KSPLogResidualHistory(ksp,dp);
 95:       KSPMonitor(ksp,i+1,0.0);
 96:       break;
 97:     }
 98:     omega = d1 / d2;                               /*   w <- (t's) / (t't) */
 99:     ierr  = VecAXPY(&alpha,P,X);     /*   x <- x + a p       */
100:     ierr  = VecAXPY(&omega,S,X);     /*   x <- x + w s       */
101:     tmp   = -omega;
102:     ierr  = VecWAXPY(&tmp,T,S,R);    /*   r <- s - w t       */
103:     if (ksp->normtype != KSP_NO_NORM) {
104:       VecNorm(R,NORM_2,&dp);
105:     }

107:     rhoold   = rho;
108:     omegaold = omega;

110:     PetscObjectTakeAccess(ksp);
111:     ksp->its++;
112:     ksp->rnorm = dp;
113:     PetscObjectGrantAccess(ksp);
114:     KSPLogResidualHistory(ksp,dp);
115:     KSPMonitor(ksp,i+1,dp);
116:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
117:     if (ksp->reason) break;
118:   }
119:   if (i == maxit) {
120:     ksp->reason = KSP_DIVERGED_ITS;
121:     i--;
122:   }
123:   *its = i+1;

125:   KSPUnwindPreconditioner(ksp,X,T);
126:   return(0);
127: }

129: EXTERN_C_BEGIN
130: int KSPCreate_BCGS(KSP ksp)
131: {
133:   ksp->data                 = (void*)0;
134:   ksp->pc_side              = PC_LEFT;
135:   ksp->ops->setup           = KSPSetUp_BCGS;
136:   ksp->ops->solve           = KSPSolve_BCGS;
137:   ksp->ops->destroy         = KSPDefaultDestroy;
138:   ksp->ops->buildsolution   = KSPDefaultBuildSolution;
139:   ksp->ops->buildresidual   = KSPDefaultBuildResidual;
140:   ksp->ops->setfromoptions  = 0;
141:   ksp->ops->view            = 0;
142:   return(0);
143: }
144: EXTERN_C_END