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

 15: static int KSPSetUp_BCGS(KSP ksp)
 16: {

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

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


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

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

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

 64:   /* Make the initial Rp == R */
 65:   VecCopy(R,RP);

 67:   rhoold   = 1.0;
 68:   alpha    = 1.0;
 69:   omegaold = 1.0;
 70:   VecSet(&zero,P);
 71:   VecSet(&zero,V);

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

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

111:     rhoold   = rho;
112:     omegaold = omega;

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

128:   KSPUnwindPreconditioner(ksp,X,T);
129:   return(0);
130: }

132: EXTERN_C_BEGIN
135: int KSPCreate_BCGS(KSP ksp)
136: {
138:   ksp->data                 = (void*)0;
139:   ksp->pc_side              = PC_LEFT;
140:   ksp->ops->setup           = KSPSetUp_BCGS;
141:   ksp->ops->solve           = KSPSolve_BCGS;
142:   ksp->ops->destroy         = KSPDefaultDestroy;
143:   ksp->ops->buildsolution   = KSPDefaultBuildSolution;
144:   ksp->ops->buildresidual   = KSPDefaultBuildResidual;
145:   ksp->ops->setfromoptions  = 0;
146:   ksp->ops->view            = 0;
147:   return(0);
148: }
149: EXTERN_C_END