Actual source code: bcgs.c

 2:  #include src/ksp/ksp/kspimpl.h

  6: static PetscErrorCode KSPSetUp_BCGS(KSP ksp)
  7: {

 11:   if (ksp->pc_side == PC_SYMMETRIC) {
 12:     SETERRQ(PETSC_ERR_SUP,"no symmetric preconditioning for KSPBCGS");
 13:   }
 14:   KSPDefaultGetWork(ksp,6);
 15:   return(0);
 16: }

 20: static PetscErrorCode  KSPSolve_BCGS(KSP ksp)
 21: {
 23:   PetscInt       i;
 24:   PetscScalar    rho,rhoold,alpha,beta,omega,omegaold,d1,d2,zero = 0.0,tmp;
 25:   Vec            X,B,V,P,R,RP,T,S;
 26:   PetscReal      dp = 0.0;


 30:   X       = ksp->vec_sol;
 31:   B       = ksp->vec_rhs;
 32:   R       = ksp->work[0];
 33:   RP      = ksp->work[1];
 34:   V       = ksp->work[2];
 35:   T       = ksp->work[3];
 36:   S       = ksp->work[4];
 37:   P       = ksp->work[5];

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

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

 55:   /* Make the initial Rp == R */
 56:   VecCopy(R,RP);

 58:   rhoold   = 1.0;
 59:   alpha    = 1.0;
 60:   omegaold = 1.0;
 61:   VecSet(&zero,P);
 62:   VecSet(&zero,V);

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

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

111:     PetscObjectTakeAccess(ksp);
112:     ksp->its++;
113:     ksp->rnorm = dp;
114:     PetscObjectGrantAccess(ksp);
115:     KSPLogResidualHistory(ksp,dp);
116:     KSPMonitor(ksp,i+1,dp);
117:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
118:     if (ksp->reason) break;
119:     i++;
120:   } while (i<ksp->max_it);

122:   if (i >= ksp->max_it) {
123:     ksp->reason = KSP_DIVERGED_ITS;
124:   }

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

130: /*MC
131:      KSPBCGS - Implements the BiCGStab (Stabilized version of BiConjugate Gradient Squared) method.

133:    Options Database Keys:
134: .   see KSPSolve()

136:    Level: beginner

138:    Notes: Reference: van der Vorst, SIAM J. Sci. Stat. Comput., 1992.

140: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG
141: M*/
145: PetscErrorCode KSPCreate_BCGS(KSP ksp)
146: {
148:   ksp->data                 = (void*)0;
149:   ksp->pc_side              = PC_LEFT;
150:   ksp->ops->setup           = KSPSetUp_BCGS;
151:   ksp->ops->solve           = KSPSolve_BCGS;
152:   ksp->ops->destroy         = KSPDefaultDestroy;
153:   ksp->ops->buildsolution   = KSPDefaultBuildSolution;
154:   ksp->ops->buildresidual   = KSPDefaultBuildResidual;
155:   ksp->ops->setfromoptions  = 0;
156:   ksp->ops->view            = 0;
157:   return(0);
158: }