Actual source code: bcgs.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>       /*I  "petscksp.h"  I*/

  6: PetscErrorCode KSPSetFromOptions_BCGS(PetscOptionItems *PetscOptionsObject,KSP ksp)
  7: {

 11:   PetscOptionsHead(PetscOptionsObject,"KSP BCGS Options");
 12:   PetscOptionsTail();
 13:   return(0);
 14: }

 18: PetscErrorCode KSPSetUp_BCGS(KSP ksp)
 19: {

 23:   KSPSetWorkVecs(ksp,6);
 24:   return(0);
 25: }


 30: PetscErrorCode KSPSolve_BCGS(KSP ksp)
 31: {
 33:   PetscInt       i;
 34:   PetscScalar    rho,rhoold,alpha,beta,omega,omegaold,d1;
 35:   Vec            X,B,V,P,R,RP,T,S;
 36:   PetscReal      dp    = 0.0,d2;
 37:   KSP_BCGS       *bcgs = (KSP_BCGS*)ksp->data;

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

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

 52:   /* with right preconditioning need to save initial guess to add to final solution */
 53:   if (ksp->pc_side == PC_RIGHT && !ksp->guess_zero) {
 54:     if (!bcgs->guess) {
 55:       VecDuplicate(X,&bcgs->guess);
 56:     }
 57:     VecCopy(X,bcgs->guess);
 58:     VecSet(X,0.0);
 59:   }

 61:   /* Test for nothing to do */
 62:   if (ksp->normtype != KSP_NORM_NONE) {
 63:     VecNorm(R,NORM_2,&dp);
 64:   }
 65:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 66:   ksp->its   = 0;
 67:   ksp->rnorm = dp;
 68:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 69:   KSPLogResidualHistory(ksp,dp);
 70:   KSPMonitor(ksp,0,dp);
 71:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 72:   if (ksp->reason) {
 73:     if (bcgs->guess) {
 74:       VecAXPY(X,1.0,bcgs->guess);
 75:     }
 76:     return(0);
 77:   }

 79:   /* Make the initial Rp == R */
 80:   VecCopy(R,RP);

 82:   rhoold   = 1.0;
 83:   alpha    = 1.0;
 84:   omegaold = 1.0;
 85:   VecSet(P,0.0);
 86:   VecSet(V,0.0);

 88:   i=0;
 89:   do {
 90:     VecDot(R,RP,&rho);       /*   rho <- (r,rp)      */
 91:     beta = (rho/rhoold) * (alpha/omegaold);
 92:     VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V);  /* p <- r - omega * beta* v + beta * p */
 93:     KSP_PCApplyBAorAB(ksp,P,V,T);  /*   v <- K p           */
 94:     VecDot(V,RP,&d1);
 95:     if (d1 == 0.0) {
 96:       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged due to Nan or Inf inner product");
 97:       else {
 98:         ksp->reason = KSP_DIVERGED_NANORINF;
 99:         break;
100:       }
101:     }
102:     alpha = rho / d1;                 /*   a <- rho / (v,rp)  */
103:     VecWAXPY(S,-alpha,V,R);     /*   s <- r - a v       */
104:     KSP_PCApplyBAorAB(ksp,S,T,R); /*   t <- K s    */
105:     VecDotNorm2(S,T,&d1,&d2);
106:     if (d2 == 0.0) {
107:       /* t is 0.  if s is 0, then alpha v == r, and hence alpha p
108:          may be our solution.  Give it a try? */
109:       VecDot(S,S,&d1);
110:       if (d1 != 0.0) {
111:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
112:         break;
113:       }
114:       VecAXPY(X,alpha,P);   /*   x <- x + a p       */
115:       PetscObjectSAWsTakeAccess((PetscObject)ksp);
116:       ksp->its++;
117:       ksp->rnorm  = 0.0;
118:       ksp->reason = KSP_CONVERGED_RTOL;
119:       PetscObjectSAWsGrantAccess((PetscObject)ksp);
120:       KSPLogResidualHistory(ksp,dp);
121:       KSPMonitor(ksp,i+1,0.0);
122:       break;
123:     }
124:     omega = d1 / d2;                               /*   w <- (t's) / (t't) */
125:     VecAXPBYPCZ(X,alpha,omega,1.0,P,S); /* x <- alpha * p + omega * s + x */
126:     VecWAXPY(R,-omega,T,S);     /*   r <- s - w t       */
127:     if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i+2) {
128:       VecNorm(R,NORM_2,&dp);
129:     }

131:     rhoold   = rho;
132:     omegaold = omega;

134:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
135:     ksp->its++;
136:     ksp->rnorm = dp;
137:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
138:     KSPLogResidualHistory(ksp,dp);
139:     KSPMonitor(ksp,i+1,dp);
140:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
141:     if (ksp->reason) break;
142:     if (rho == 0.0) {
143:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
144:       break;
145:     }
146:     i++;
147:   } while (i<ksp->max_it);

149:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;

151:   KSPUnwindPreconditioner(ksp,X,T);
152:   if (bcgs->guess) {
153:     VecAXPY(X,1.0,bcgs->guess);
154:   }
155:   return(0);
156: }

160: PetscErrorCode KSPBuildSolution_BCGS(KSP ksp,Vec v,Vec *V)
161: {
163:   KSP_BCGS       *bcgs = (KSP_BCGS*)ksp->data;

166:   if (ksp->pc_side == PC_RIGHT) {
167:     if (v) {
168:       KSP_PCApply(ksp,ksp->vec_sol,v);
169:       if (bcgs->guess) {
170:         VecAXPY(v,1.0,bcgs->guess);
171:       }
172:       *V = v;
173:     } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
174:   } else {
175:     if (v) {
176:       VecCopy(ksp->vec_sol,v); *V = v;
177:     } else *V = ksp->vec_sol;
178:   }
179:   return(0);
180: }

184: PetscErrorCode KSPReset_BCGS(KSP ksp)
185: {
186:   KSP_BCGS       *cg = (KSP_BCGS*)ksp->data;

190:   VecDestroy(&cg->guess);
191:   return(0);
192: }

196: PetscErrorCode KSPDestroy_BCGS(KSP ksp)
197: {

201:   KSPReset_BCGS(ksp);
202:   KSPDestroyDefault(ksp);
203:   return(0);
204: }

206: /*MC
207:      KSPBCGS - Implements the BiCGStab (Stabilized version of BiConjugate Gradient) method.

209:    Options Database Keys:
210: .   see KSPSolve()

212:    Level: beginner

214:    Notes: See KSPBCGSL for additional stabilization
215:           Supports left and right preconditioning but not symmetric

217:    References:
218: .    1. -   van der Vorst, SIAM J. Sci. Stat. Comput., 1992.

220: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL, KSPFBICG, KSPSetPCSide()
221: M*/
224: PETSC_EXTERN PetscErrorCode KSPCreate_BCGS(KSP ksp)
225: {
227:   KSP_BCGS       *bcgs;

230:   PetscNewLog(ksp,&bcgs);

232:   ksp->data                = bcgs;
233:   ksp->ops->setup          = KSPSetUp_BCGS;
234:   ksp->ops->solve          = KSPSolve_BCGS;
235:   ksp->ops->destroy        = KSPDestroy_BCGS;
236:   ksp->ops->reset          = KSPReset_BCGS;
237:   ksp->ops->buildsolution  = KSPBuildSolution_BCGS;
238:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
239:   ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;

241:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
242:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
243:   return(0);
244: }