Actual source code: fbcgs.c

petsc-3.3-p5 2012-12-01
  2: /*
  3:     This file implements flexible BiCGStab contributed by Jie Chen.
  4:     Only right preconditioning is supported. 
  5: */
  6: #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>       /*I  "petscksp.h"  I*/

  8: /* copied from KSPBuildSolution_GCR() */
 11: PetscErrorCode  KSPBuildSolution_FBCGS(KSP ksp, Vec v, Vec *V)
 12: {
 14:   Vec            x;
 15: 
 17:   x = ksp->vec_sol;
 18:   if (v) {
 19:     VecCopy(x, v);
 20:     if (V) *V = v;
 21:   } else if (V) {
 22:     *V = ksp->vec_sol;
 23:   }
 24:   return(0);
 25: }

 29: static PetscErrorCode KSPSetUp_FBCGS(KSP ksp)
 30: {
 32: 
 34:   KSPDefaultGetWork(ksp,8);
 35:   return(0);
 36: }

 38: /* Only need a few hacks from KSPSolve_BCGS */
 39: #include <petsc-private/pcimpl.h>            /*I "petscksp.h" I*/
 42: static PetscErrorCode  KSPSolve_FBCGS(KSP ksp)
 43: {
 45:   PetscInt       i;
 46:   PetscScalar    rho,rhoold,alpha,beta,omega,omegaold,d1,d2;
 47:   Vec            X,B,V,P,R,RP,T,S,P2,S2;
 48:   PetscReal      dp = 0.0;
 49:   KSP_BCGS       *bcgs = (KSP_BCGS*)ksp->data;
 50:   PC             pc;

 53:   X       = ksp->vec_sol;
 54:   B       = ksp->vec_rhs;
 55:   R       = ksp->work[0];
 56:   RP      = ksp->work[1];
 57:   V       = ksp->work[2];
 58:   T       = ksp->work[3];
 59:   S       = ksp->work[4];
 60:   P       = ksp->work[5];
 61:   S2      = ksp->work[6];
 62:   P2      = ksp->work[7];

 64:   /* Only supports right preconditioning */
 65:   if (ksp->pc_side != PC_RIGHT)
 66:     SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"KSP fbcgs does not support %s",PCSides[ksp->pc_side]);
 67:   if (!ksp->guess_zero) {
 68:     if (!bcgs->guess) {
 69:       VecDuplicate(X,&bcgs->guess);
 70:     }
 71:     VecCopy(X,bcgs->guess);
 72:   } else {
 73:     VecSet(X,0.0);
 74:   }

 76:   /* Compute initial residual */
 77:   KSPGetPC(ksp,&pc);
 78:   if (pc->setupcalled < 2) {
 79:     PCSetUp(pc);
 80:   }
 81:   if (!ksp->guess_zero) {
 82:     MatMult(pc->mat,X,S2);
 83:     VecCopy(B,R);
 84:     VecAXPY(R,-1.0,S2);
 85:   } else {
 86:     VecCopy(B,R);
 87:   }

 89:   /* Test for nothing to do */
 90:   if (ksp->normtype != KSP_NORM_NONE) {
 91:     VecNorm(R,NORM_2,&dp);
 92:   }
 93:   PetscObjectTakeAccess(ksp);
 94:   ksp->its   = 0;
 95:   ksp->rnorm = dp;
 96:   PetscObjectGrantAccess(ksp);
 97:   KSPLogResidualHistory(ksp,dp);
 98:   KSPMonitor(ksp,0,dp);
 99:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
100:   if (ksp->reason) return(0);

102:   /* Make the initial Rp == R */
103:   VecCopy(R,RP);

105:   rhoold   = 1.0;
106:   alpha    = 1.0;
107:   omegaold = 1.0;
108:   VecSet(P,0.0);
109:   VecSet(V,0.0);
110: 
111:   i=0;
112:   do {
113:     VecDot(R,RP,&rho); /* rho <- (r,rp) */
114:     beta = (rho/rhoold) * (alpha/omegaold);
115:     VecAXPBYPCZ(P,1.0,-omegaold*beta,beta,R,V); /* p <- r - omega * beta* v + beta * p */

117:     if (pc->setupcalled < 2) {
118:       PCSetUp(pc);
119:     }
120:     PCApply(pc,P,P2); /* p2 <- K p */
121:     MatMult(pc->mat,P2,V); /* v <- A p2 */
122: 
123:     VecDot(V,RP,&d1);
124:     if (d1 == 0.0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Divide by zero");
125:     alpha = rho / d1; /* alpha <- rho / (v,rp) */
126:     VecWAXPY(S,-alpha,V,R);  /* s <- r - alpha v */

128:     if (pc->setupcalled < 2) {
129:       PCSetUp(pc);
130:     }
131:     PCApply(pc,S,S2); /* s2 <- K s */
132:     MatMult(pc->mat,S2,T); /* t <- A s2 */

134:     VecDotNorm2(S,T,&d1,&d2);
135:     if (d2 == 0.0) {
136:       /* t is 0. if s is 0, then alpha v == r, and hence alpha p may be our solution. Give it a try? */
137:       VecDot(S,S,&d1);
138:       if (d1 != 0.0) {
139:         ksp->reason = KSP_DIVERGED_BREAKDOWN;
140:         break;
141:       }
142:       VecAXPY(X,alpha,P2);   /* x <- x + alpha p2 */
143:       PetscObjectTakeAccess(ksp);
144:       ksp->its++;
145:       ksp->rnorm  = 0.0;
146:       ksp->reason = KSP_CONVERGED_RTOL;
147:       PetscObjectGrantAccess(ksp);
148:       KSPLogResidualHistory(ksp,dp);
149:       KSPMonitor(ksp,i+1,0.0);
150:       break;
151:     }
152:     omega = d1 / d2; /* omega <- (t's) / (t't) */
153:     VecAXPBYPCZ(X,alpha,omega,1.0,P2,S2); /* x <- alpha * p2 + omega * s2 + x */

155:     VecWAXPY(R,-omega,T,S); /* r <- s - omega t */
156:     if (ksp->normtype != KSP_NORM_NONE && ksp->chknorm < i+2) {
157:       VecNorm(R,NORM_2,&dp);
158:     }

160:     rhoold   = rho;
161:     omegaold = omega;

163:     PetscObjectTakeAccess(ksp);
164:     ksp->its++;
165:     ksp->rnorm = dp;
166:     PetscObjectGrantAccess(ksp);
167:     KSPLogResidualHistory(ksp,dp);
168:     KSPMonitor(ksp,i+1,dp);
169:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
170:     if (ksp->reason) break;
171:     if (rho == 0.0) {
172:       ksp->reason = KSP_DIVERGED_BREAKDOWN;
173:       break;
174:     }
175:     i++;
176:   } while (i<ksp->max_it);

178:   if (i >= ksp->max_it) {
179:     ksp->reason = KSP_DIVERGED_ITS;
180:   }
181:   return(0);
182: }

184: /*MC
185:      KSPFBCGS - Implements flexible BiCGStab (Stabilized version of BiConjugate Gradient Squared) method.

187:    Options Database Keys:
188: .   see KSPSolve()

190:    Level: beginner

192:    Notes: Only supports right preconditioning 

194: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPFBCGSL, KSPSetPCSide()
195: M*/
196: EXTERN_C_BEGIN
199: PetscErrorCode  KSPCreate_FBCGS(KSP ksp)
200: {
202:   KSP_BCGS       *bcgs;

205:   PetscNewLog(ksp,KSP_BCGS,&bcgs);
206:   ksp->data                 = bcgs;
207:   ksp->ops->setup           = KSPSetUp_FBCGS;
208:   ksp->ops->solve           = KSPSolve_FBCGS;
209:   ksp->ops->destroy         = KSPDestroy_BCGS;
210:   ksp->ops->reset           = KSPReset_BCGS;
211:   ksp->ops->buildsolution   = KSPBuildSolution_FBCGS;
212:   ksp->ops->buildresidual   = KSPDefaultBuildResidual;
213:   ksp->ops->setfromoptions  = KSPSetFromOptions_BCGS;
214:   ksp->ops->view            = KSPView_BCGS;
215:   ksp->pc_side              = PC_RIGHT; /* set default PC side */
216: 
217:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
218:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
219:   return(0);
220: }
221: EXTERN_C_END