Actual source code: bicg.c

petsc-master 2020-07-13
Report Typos and Errors

  2:  #include <petsc/private/kspimpl.h>

  4: static PetscErrorCode KSPSetUp_BiCG(KSP ksp)
  5: {

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

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

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

 30:   X  = ksp->vec_sol;
 31:   B  = ksp->vec_rhs;
 32:   Rl = ksp->work[0];
 33:   Zl = ksp->work[1];
 34:   Pl = ksp->work[2];
 35:   Rr = ksp->work[3];
 36:   Zr = ksp->work[4];
 37:   Pr = ksp->work[5];

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

 41:   if (!ksp->guess_zero) {
 42:     KSP_MatMult(ksp,Amat,X,Rr);      /*   r <- b - Ax       */
 43:     VecAYPX(Rr,-1.0,B);
 44:   } else {
 45:     VecCopy(B,Rr);           /*     r <- b (x is 0) */
 46:   }
 47:   VecCopy(Rr,Rl);
 48:   KSP_PCApply(ksp,Rr,Zr);     /*     z <- Br         */
 49:   VecConjugate(Rl);
 50:   KSP_PCApplyTranspose(ksp,Rl,Zl);
 51:   VecConjugate(Rl);
 52:   VecConjugate(Zl);
 53:   if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
 54:     VecNorm(Zr,NORM_2,&dp);  /*    dp <- z'*z       */
 55:   } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
 56:     VecNorm(Rr,NORM_2,&dp);  /*    dp <- r'*r       */
 57:   } else dp = 0.0;

 59:   KSPCheckNorm(ksp,dp);
 60:   KSPMonitor(ksp,0,dp);
 61:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 62:   ksp->its   = 0;
 63:   ksp->rnorm = dp;
 64:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 65:   KSPLogResidualHistory(ksp,dp);
 66:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 67:   if (ksp->reason) return(0);

 69:   i = 0;
 70:   do {
 71:     VecDot(Zr,Rl,&beta);       /*     beta <- r'z     */
 72:     KSPCheckDot(ksp,beta);
 73:     if (!i) {
 74:       if (beta == 0.0) {
 75:         ksp->reason = KSP_DIVERGED_BREAKDOWN_BICG;
 76:         return(0);
 77:       }
 78:       VecCopy(Zr,Pr);       /*     p <- z          */
 79:       VecCopy(Zl,Pl);
 80:     } else {
 81:       b    = beta/betaold;
 82:       VecAYPX(Pr,b,Zr);  /*     p <- z + b* p   */
 83:       b    = PetscConj(b);
 84:       VecAYPX(Pl,b,Zl);
 85:     }
 86:     betaold = beta;
 87:     KSP_MatMult(ksp,Amat,Pr,Zr); /*     z <- Kp         */
 88:     VecConjugate(Pl);
 89:     KSP_MatMultTranspose(ksp,Amat,Pl,Zl);
 90:     VecConjugate(Pl);
 91:     VecConjugate(Zl);
 92:     VecDot(Zr,Pl,&dpi);            /*     dpi <- z'p      */
 93:     KSPCheckDot(ksp,dpi);
 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 if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
108:       VecNorm(Rr,NORM_2,&dp);  /*    dp <- r'*r       */
109:     } else dp = 0.0;

111:     KSPCheckNorm(ksp,dp);
112:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
113:     ksp->its   = i+1;
114:     ksp->rnorm = dp;
115:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
116:     KSPLogResidualHistory(ksp,dp);
117:     KSPMonitor(ksp,i+1,dp);
118:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
119:     if (ksp->reason) break;
120:     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
121:       KSP_PCApply(ksp,Rr,Zr);  /* z <- Br  */
122:       VecConjugate(Rl);
123:       KSP_PCApplyTranspose(ksp,Rl,Zl);
124:       VecConjugate(Rl);
125:       VecConjugate(Zl);
126:     }
127:     i++;
128:   } while (i<ksp->max_it);
129:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
130:   return(0);
131: }

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

137:    Options Database Keys:
138: .   see KSPSolve()

140:    Level: beginner

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

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

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

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

158:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
159:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
160:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

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