Actual source code: bicg.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petsc-private/kspimpl.h>

  6: PetscErrorCode KSPSetUp_BiCG(KSP ksp)
  7: {

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

 20: PetscErrorCode  KSPSolve_BiCG(KSP ksp)
 21: {
 23:   PetscInt       i;
 24:   PetscBool      diagonalscale;
 25:   PetscScalar    dpi,a=1.0,beta,betaold=1.0,b,ma;
 26:   PetscReal      dp;
 27:   Vec            X,B,Zl,Zr,Rl,Rr,Pl,Pr;
 28:   Mat            Amat,Pmat;

 31:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 32:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 34:   X  = ksp->vec_sol;
 35:   B  = ksp->vec_rhs;
 36:   Rl = ksp->work[0];
 37:   Zl = ksp->work[1];
 38:   Pl = ksp->work[2];
 39:   Rr = ksp->work[3];
 40:   Zr = ksp->work[4];
 41:   Pr = ksp->work[5];

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

 45:   if (!ksp->guess_zero) {
 46:     KSP_MatMult(ksp,Amat,X,Rr);      /*   r <- b - Ax       */
 47:     VecAYPX(Rr,-1.0,B);
 48:   } else {
 49:     VecCopy(B,Rr);           /*     r <- b (x is 0) */
 50:   }
 51:   VecCopy(Rr,Rl);
 52:   KSP_PCApply(ksp,Rr,Zr);     /*     z <- Br         */
 53:   VecConjugate(Rl);
 54:   KSP_PCApplyTranspose(ksp,Rl,Zl);
 55:   VecConjugate(Rl);
 56:   VecConjugate(Zl);
 57:   if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
 58:     VecNorm(Zr,NORM_2,&dp);  /*    dp <- z'*z       */
 59:   } else {
 60:     VecNorm(Rr,NORM_2,&dp);  /*    dp <- r'*r       */
 61:   }
 62:   KSPMonitor(ksp,0,dp);
 63:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 64:   ksp->its   = 0;
 65:   ksp->rnorm = dp;
 66:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 67:   KSPLogResidualHistory(ksp,dp);
 68:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 69:   if (ksp->reason) return(0);

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

131: /*MC
132:      KSPBICG - Implements the Biconjugate gradient method (similar to running the conjugate
133:          gradient on the normal equations).

135:    Options Database Keys:
136: .   see KSPSolve()

138:    Level: beginner

140:    Notes: this method requires that one be apply to apply the transpose of the preconditioner and operator
141:          as well as the operator and preconditioner.
142:          Supports only left preconditioning

144:          See KSPCGNE for code that EXACTLY runs the preconditioned conjugate gradient method on the
145:          normal equations

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

149: M*/
152: PETSC_EXTERN PetscErrorCode KSPCreate_BiCG(KSP ksp)
153: {

157:   ksp->data = (void*)0;
158:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
159:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);

161:   ksp->ops->setup          = KSPSetUp_BiCG;
162:   ksp->ops->solve          = KSPSolve_BiCG;
163:   ksp->ops->destroy        = KSPDestroyDefault;
164:   ksp->ops->view           = 0;
165:   ksp->ops->setfromoptions = 0;
166:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
167:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
168:   return(0);
169: }