Actual source code: bcgs.c

  2: #include <private/kspimpl.h>

  4: typedef struct {
  5:   Vec guess;        /* if using right preconditioning with nonzero initial guess must keep that around to "fix" solution */
  6: } KSP_BCGS;

 10: static PetscErrorCode KSPSetUp_BCGS(KSP ksp)
 11: {

 15:   KSPDefaultGetWork(ksp,6);
 16:   return(0);
 17: }

 21: static PetscErrorCode  KSPSolve_BCGS(KSP ksp)
 22: {
 24:   PetscInt       i;
 25:   PetscScalar    rho,rhoold,alpha,beta,omega,omegaold,d1,d2;
 26:   Vec            X,B,V,P,R,RP,T,S;
 27:   PetscReal      dp = 0.0;
 28:   KSP_BCGS       *bcgs = (KSP_BCGS*)ksp->data;

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

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

 43:   /* with right preconditioning need to save initial guess to add to final solution */
 44:   if (ksp->pc_side == PC_RIGHT && !ksp->guess_zero) {
 45:     if (!bcgs->guess) {
 46:       VecDuplicate(X,&bcgs->guess);
 47:     }
 48:     VecCopy(X,bcgs->guess);
 49:     VecSet(X,0.0);
 50:   }

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

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

 68:   rhoold   = 1.0;
 69:   alpha    = 1.0;
 70:   omegaold = 1.0;
 71:   VecSet(P,0.0);
 72:   VecSet(V,0.0);

 74:   i=0;
 75:   do {
 76:     VecDot(R,RP,&rho);       /*   rho <- (r,rp)      */
 77:     beta = (rho/rhoold) * (alpha/omegaold);
 78:     VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V);  /* p <- r - omega * beta* v + beta * p */
 79:     KSP_PCApplyBAorAB(ksp,P,V,T);  /*   v <- K p           */
 80:     VecDot(V,RP,&d1);
 81:     if (d1 == 0.0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Divide by zero");
 82:     alpha = rho / d1;                 /*   a <- rho / (v,rp)  */
 83:     VecWAXPY(S,-alpha,V,R);      /*   s <- r - a v       */
 84:     KSP_PCApplyBAorAB(ksp,S,T,R);/*   t <- K s    */
 85:     VecDotNorm2(S,T,&d1,&d2);
 86:     if (d2 == 0.0) {
 87:       /* t is 0.  if s is 0, then alpha v == r, and hence alpha p
 88:          may be our solution.  Give it a try? */
 89:       VecDot(S,S,&d1);
 90:       if (d1 != 0.0) {
 91:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
 92:         break;
 93:       }
 94:       VecAXPY(X,alpha,P);   /*   x <- x + a p       */
 95:       PetscObjectTakeAccess(ksp);
 96:       ksp->its++;
 97:       ksp->rnorm  = 0.0;
 98:       ksp->reason = KSP_CONVERGED_RTOL;
 99:       PetscObjectGrantAccess(ksp);
100:       KSPLogResidualHistory(ksp,dp);
101:       KSPMonitor(ksp,i+1,0.0);
102:       break;
103:     }
104:     omega = d1 / d2;                               /*   w <- (t's) / (t't) */
105:     VecAXPBYPCZ(X,alpha,omega,1.0,P,S); /* x <- alpha * p + omega * s + x */
106:     VecWAXPY(R,-omega,T,S);     /*   r <- s - w t       */
107:     if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i+2) {
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:     if (rho == 0.0) {
123:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
124:       break;
125:     }
126:     i++;
127:   } while (i<ksp->max_it);

129:   if (i >= ksp->max_it) {
130:     ksp->reason = KSP_DIVERGED_ITS;
131:   }

133:   KSPUnwindPreconditioner(ksp,X,T);
134:   if (bcgs->guess) {
135:     VecAXPY(X,1.0,bcgs->guess);
136:   }
137:   return(0);
138: }

142: PetscErrorCode KSPBuildSolution_BCGS(KSP ksp,Vec v,Vec *V)
143: {
145:   KSP_BCGS       *bcgs = (KSP_BCGS*)ksp->data;

148:   if (ksp->pc_side == PC_RIGHT) {
149:     if (v) {
150:       KSP_PCApply(ksp,ksp->vec_sol,v);
151:       if (bcgs->guess) {
152:         VecAXPY(v,1.0,bcgs->guess);
153:       }
154:       *V = v;
155:     } else SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Not working with right preconditioner");
156:   } else {
157:     if (v) {VecCopy(ksp->vec_sol,v); *V = v;}
158:     else { *V = ksp->vec_sol; }
159:   }
160:   return(0);
161: }

165: PetscErrorCode KSPReset_BCGS(KSP ksp)
166: {
167:   KSP_BCGS       *cg = (KSP_BCGS*)ksp->data;

171:   VecDestroy(&cg->guess);
172:   return(0);
173: }

177: PetscErrorCode KSPDestroy_BCGS(KSP ksp)
178: {

182:   KSPReset_BCGS(ksp);
183:   KSPDefaultDestroy(ksp);
184:   return(0);
185: }

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

190:    Options Database Keys:
191: .   see KSPSolve()

193:    Level: beginner

195:    Notes: See KSPBCGSL for additional stabilization
196:           Supports left and right preconditioning but not symmetric

198:    References: van der Vorst, SIAM J. Sci. Stat. Comput., 1992.


201: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL, KSPSetPCSide()
202: M*/
206: PetscErrorCode  KSPCreate_BCGS(KSP ksp)
207: {

211:   PetscNewLog(ksp,KSP_BCGS,&ksp->data);
212:   ksp->ops->setup           = KSPSetUp_BCGS;
213:   ksp->ops->solve           = KSPSolve_BCGS;
214:   ksp->ops->destroy         = KSPDestroy_BCGS;
215:   ksp->ops->reset           = KSPReset_BCGS;
216:   ksp->ops->buildsolution   = KSPBuildSolution_BCGS;
217:   ksp->ops->buildresidual   = KSPDefaultBuildResidual;
218:   ksp->ops->setfromoptions  = 0;
219:   ksp->ops->view            = 0;

221:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
222:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
223:   return(0);
224: }