Actual source code: bicg.c

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

  6: PetscErrorCode KSPSetUp_BiCG(KSP ksp)
  7: {

 11:   /* check user parameters and functions */
 12:   if (ksp->pc_side == PC_RIGHT) {
 13:     SETERRQ(2,"no right preconditioning for KSPBiCG");
 14:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 15:     SETERRQ(2,"no symmetric preconditioning for KSPBiCG");
 16:   }

 18:   /* get work vectors from user code */
 19:   KSPDefaultGetWork(ksp,6);
 20:   return(0);
 21: }

 25: PetscErrorCode  KSPSolve_BiCG(KSP ksp)
 26: {
 28:   PetscInt       i;
 29:   PetscTruth     diagonalscale;
 30:   PetscScalar    dpi,a=1.0,beta,betaold=1.0,b,mone=-1.0,ma;
 31:   PetscReal      dp;
 32:   Vec            X,B,Zl,Zr,Rl,Rr,Pl,Pr;
 33:   Mat            Amat,Pmat;
 34:   MatStructure   pflag;

 37:   PCDiagonalScale(ksp->pc,&diagonalscale);
 38:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);

 40:   X       = ksp->vec_sol;
 41:   B       = ksp->vec_rhs;
 42:   Rl      = ksp->work[0];
 43:   Zl      = ksp->work[1];
 44:   Pl      = ksp->work[2];
 45:   Rr      = ksp->work[3];
 46:   Zr      = ksp->work[4];
 47:   Pr      = ksp->work[5];

 49:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);

 51:   if (!ksp->guess_zero) {
 52:     KSP_MatMult(ksp,Amat,X,Rr);      /*   r <- b - Ax       */
 53:     VecAYPX(&mone,B,Rr);
 54:   } else {
 55:     VecCopy(B,Rr);           /*     r <- b (x is 0) */
 56:   }
 57:   VecCopy(Rr,Rl);
 58:   KSP_PCApply(ksp,Rr,Zr);     /*     z <- Br         */
 59:   VecConjugate(Rl);
 60:   KSP_PCApplyTranspose(ksp,Rl,Zl);
 61:   VecConjugate(Rl);
 62:   VecConjugate(Zl);
 63:   if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
 64:     VecNorm(Zr,NORM_2,&dp);  /*    dp <- z'*z       */
 65:   } else {
 66:     VecNorm(Rr,NORM_2,&dp);  /*    dp <- r'*r       */
 67:   }
 68:   KSPMonitor(ksp,0,dp);
 69:   PetscObjectTakeAccess(ksp);
 70:   ksp->its   = 0;
 71:   ksp->rnorm = dp;
 72:   PetscObjectGrantAccess(ksp);
 73:   KSPLogResidualHistory(ksp,dp);
 74:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 75:   if (ksp->reason) return(0);

 77:   i = 0;
 78:   do {
 79:      VecDot(Zr,Rl,&beta);       /*     beta <- r'z     */
 80:      if (!i) {
 81:        if (beta == 0.0) {
 82:          ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
 83:          return(0);
 84:        }
 85:        VecCopy(Zr,Pr);       /*     p <- z          */
 86:        VecCopy(Zl,Pl);
 87:      } else {
 88:        b = beta/betaold;
 89:        VecAYPX(&b,Zr,Pr);  /*     p <- z + b* p   */
 90:        b = PetscConj(b);
 91:        VecAYPX(&b,Zl,Pl);
 92:      }
 93:      betaold = beta;
 94:      KSP_MatMult(ksp,Amat,Pr,Zr);    /*     z <- Kp         */
 95:      VecConjugate(Pl);
 96:      KSP_MatMultTranspose(ksp,Amat,Pl,Zl);
 97:      VecConjugate(Pl);
 98:      VecConjugate(Zl);
 99:      VecDot(Zr,Pl,&dpi);               /*     dpi <- z'p      */
100:      a = beta/dpi;                                 /*     a = beta/p'z    */
101:      VecAXPY(&a,Pr,X);       /*     x <- x + ap     */
102:      ma = -a;
103:      VecAXPY(&ma,Zr,Rr);CHKERRQ(ierr)
104:      ma = PetscConj(ma);
105:      VecAXPY(&ma,Zl,Rl);
106:      if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
107:        KSP_PCApply(ksp,Rr,Zr);  /*     z <- Br         */
108:        VecConjugate(Rl);
109:        KSP_PCApplyTranspose(ksp,Rl,Zl);
110:        VecConjugate(Rl);
111:        VecConjugate(Zl);
112:        VecNorm(Zr,NORM_2,&dp);  /*    dp <- z'*z       */
113:      } else {
114:        VecNorm(Rr,NORM_2,&dp);  /*    dp <- r'*r       */
115:      }
116:      PetscObjectTakeAccess(ksp);
117:      ksp->its   = i+1;
118:      ksp->rnorm = dp;
119:      PetscObjectGrantAccess(ksp);
120:      KSPLogResidualHistory(ksp,dp);
121:      KSPMonitor(ksp,i+1,dp);
122:      (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
123:      if (ksp->reason) break;
124:      if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
125:        KSP_PCApply(ksp,Rr,Zr);  /* z <- Br  */
126:        VecConjugate(Rl);
127:        KSP_PCApplyTranspose(ksp,Rl,Zl);
128:        VecConjugate(Rl);
129:        VecConjugate(Zl);
130:      }
131:      i++;
132:   } while (i<ksp->max_it);
133:   if (i >= ksp->max_it) {
134:     ksp->reason = KSP_DIVERGED_ITS;
135:   }
136:   return(0);
137: }

141: PetscErrorCode KSPDestroy_BiCG(KSP ksp)
142: {

146:   KSPDefaultFreeWork(ksp);
147:   return(0);
148: }

150: /*MC
151:      KSPBICG - Implements the Biconjugate gradient method (essentially running the conjugate
152:          gradient on the normal equations).

154:    Options Database Keys:
155: .   see KSPSolve()

157:    Level: beginner

159: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBCGS

161: M*/
165: PetscErrorCode KSPCreate_BiCG(KSP ksp)
166: {
168:   ksp->data                      = (void*)0;
169:   ksp->pc_side                   = PC_LEFT;
170:   ksp->ops->setup                = KSPSetUp_BiCG;
171:   ksp->ops->solve                = KSPSolve_BiCG;
172:   ksp->ops->destroy              = KSPDestroy_BiCG;
173:   ksp->ops->view                 = 0;
174:   ksp->ops->setfromoptions       = 0;
175:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
176:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

178:   return(0);
179: }