Actual source code: bcgs.c

petsc-3.15.0 2021-04-05
Report Typos and Errors

  2: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>

  4: PetscErrorCode KSPSetFromOptions_BCGS(PetscOptionItems *PetscOptionsObject,KSP ksp)
  5: {

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

 14: PetscErrorCode KSPSetUp_BCGS(KSP ksp)
 15: {

 19:   KSPSetWorkVecs(ksp,6);
 20:   return(0);
 21: }

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

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

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

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

 54:   /* Test for nothing to do */
 55:   if (ksp->normtype != KSP_NORM_NONE) {
 56:     VecNorm(R,NORM_2,&dp);
 57:     KSPCheckNorm(ksp,dp);
 58:   }
 59:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 60:   ksp->its   = 0;
 61:   ksp->rnorm = dp;
 62:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 63:   KSPLogResidualHistory(ksp,dp);
 64:   KSPMonitor(ksp,0,dp);
 65:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 66:   if (ksp->reason) {
 67:     if (bcgs->guess) {
 68:       VecAXPY(X,1.0,bcgs->guess);
 69:     }
 70:     return(0);
 71:   }

 73:   /* Make the initial Rp == R */
 74:   VecCopy(R,RP);

 76:   rhoold   = 1.0;
 77:   alpha    = 1.0;
 78:   omegaold = 1.0;
 79:   VecSet(P,0.0);
 80:   VecSet(V,0.0);

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

127:     rhoold   = rho;
128:     omegaold = omega;

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

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

147:   KSPUnwindPreconditioner(ksp,X,T);
148:   if (bcgs->guess) {
149:     VecAXPY(X,1.0,bcgs->guess);
150:   }
151:   return(0);
152: }

154: PetscErrorCode KSPBuildSolution_BCGS(KSP ksp,Vec v,Vec *V)
155: {
157:   KSP_BCGS       *bcgs = (KSP_BCGS*)ksp->data;

160:   if (ksp->pc_side == PC_RIGHT) {
161:     if (v) {
162:       KSP_PCApply(ksp,ksp->vec_sol,v);
163:       if (bcgs->guess) {
164:         VecAXPY(v,1.0,bcgs->guess);
165:       }
166:       *V = v;
167:     } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
168:   } else {
169:     if (v) {
170:       VecCopy(ksp->vec_sol,v); *V = v;
171:     } else *V = ksp->vec_sol;
172:   }
173:   return(0);
174: }

176: PetscErrorCode KSPReset_BCGS(KSP ksp)
177: {
178:   KSP_BCGS       *cg = (KSP_BCGS*)ksp->data;

182:   VecDestroy(&cg->guess);
183:   return(0);
184: }

186: PetscErrorCode KSPDestroy_BCGS(KSP ksp)
187: {

191:   KSPReset_BCGS(ksp);
192:   KSPDestroyDefault(ksp);
193:   return(0);
194: }

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

199:    Options Database Keys:
200: .   see KSPSolve()

202:    Level: beginner

204:    Notes:
205:     See KSPBCGSL for additional stabilization
206:           Supports left and right preconditioning but not symmetric

208:    References:
209: .    1. -   van der Vorst, SIAM J. Sci. Stat. Comput., 1992.

211: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPBCGSL, KSPFBICG, KSPSetPCSide()
212: M*/
213: PETSC_EXTERN PetscErrorCode KSPCreate_BCGS(KSP ksp)
214: {
216:   KSP_BCGS       *bcgs;

219:   PetscNewLog(ksp,&bcgs);

221:   ksp->data                = bcgs;
222:   ksp->ops->setup          = KSPSetUp_BCGS;
223:   ksp->ops->solve          = KSPSolve_BCGS;
224:   ksp->ops->destroy        = KSPDestroy_BCGS;
225:   ksp->ops->reset          = KSPReset_BCGS;
226:   ksp->ops->buildsolution  = KSPBuildSolution_BCGS;
227:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
228:   ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;

230:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
231:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
232:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
233:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
234:   return(0);
235: }