Actual source code: groppcg.c

petsc-3.4.0 2013-05-13
  1: /*
  2:  author: Pieter Ghysels, Universiteit Antwerpen, Intel Exascience lab Flanders

  4:  This file implements a CG method due to B. Gropp that can overlap one reduction
  5:  with application of the preconditioner and the other reduction can be
  6:  overlapped with the matrix-vector product.

  8:  see:
  9:  http://www.cs.uiuc.edu/~wgropp/bib/talks/tdata/2012/icerm.pdf
 10:  */

 12: #include <petsc-private/kspimpl.h>

 14: /*
 15:  KSPSetUp_GROPPCG - Sets up the workspace needed by the GROPPCG method.

 17:  This is called once, usually automatically by KSPSolve() or KSPSetUp()
 18:  but can be called directly by KSPSetUp()
 19: */
 22: PetscErrorCode KSPSetUp_GROPPCG(KSP ksp)
 23: {

 27:   KSPSetWorkVecs(ksp,6);
 28:   return(0);
 29: }

 31: /*
 32:  KSPSolve_GROPPCG

 34:  Input Parameter:
 35:  .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
 36:              example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
 37: */
 40: PetscErrorCode  KSPSolve_GROPPCG(KSP ksp)
 41: {
 43:   PetscInt       i;
 44:   PetscScalar    alpha,beta = 0.0,gamma,gammaNew,t;
 45:   PetscReal      dp = 0.0;
 46:   Vec            x,b,r,p,s,S,z,Z;
 47:   Mat            Amat,Pmat;
 48:   MatStructure   pflag;
 49:   PetscBool      diagonalscale;

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

 55:   x = ksp->vec_sol;
 56:   b = ksp->vec_rhs;
 57:   r = ksp->work[0];
 58:   p = ksp->work[1];
 59:   s = ksp->work[2];
 60:   S = ksp->work[3];
 61:   z = ksp->work[4];
 62:   Z = ksp->work[5];

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

 66:   ksp->its = 0;
 67:   if (!ksp->guess_zero) {
 68:     KSP_MatMult(ksp,Amat,x,r);            /*     r <- b - Ax     */
 69:     VecAYPX(r,-1.0,b);
 70:   } else {
 71:     VecCopy(b,r);                         /*     r <- b (x is 0) */
 72:   }

 74:   KSP_PCApply(ksp,r,z);                   /*     z <- Br   */
 75:   VecCopy(z,p);                           /*     p <- z    */
 76:   VecDotBegin(r,z,&gamma);                  /*     gamma <- z'*r       */
 77:   PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)r));
 78:   KSP_MatMult(ksp,Amat,p,s);              /*     s <- Ap   */
 79:   VecDotEnd(r,z,&gamma);                  /*     gamma <- z'*r       */

 81:   switch (ksp->normtype) {
 82:   case KSP_NORM_PRECONDITIONED:
 83:     /* This could be merged with the computation of gamma above */
 84:     VecNorm(z,NORM_2,&dp);                /*     dp <- z'*z = e'*A'*B'*B*A'*e'     */
 85:     break;
 86:   case KSP_NORM_UNPRECONDITIONED:
 87:     /* This could be merged with the computation of gamma above */
 88:     VecNorm(r,NORM_2,&dp);                /*     dp <- r'*r = e'*A'*A*e            */
 89:     break;
 90:   case KSP_NORM_NATURAL:
 91:     if (PetscIsInfOrNanScalar(gamma)) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
 92:     dp = PetscSqrtReal(PetscAbsScalar(gamma));                  /*     dp <- r'*z = r'*B*r = e'*A'*B*A*e */
 93:     break;
 94:   case KSP_NORM_NONE:
 95:     dp = 0.0;
 96:     break;
 97:   default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
 98:   }
 99:   KSPLogResidualHistory(ksp,dp);
100:   KSPMonitor(ksp,0,dp);
101:   ksp->rnorm = dp;
102:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
103:   if (ksp->reason) return(0);

105:   i = 0;
106:   do {
107:     ksp->its = i+1;
108:     i++;

110:     VecDotBegin(p,s,&t);
111:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)p));

113:     KSP_PCApply(ksp,s,S);         /*   S <- Bs       */

115:     VecDotEnd(p,s,&t);

117:     alpha = gamma / t;
118:     VecAXPY(x, alpha,p);   /*     x <- x + alpha * p   */
119:     VecAXPY(r,-alpha,s);   /*     r <- r - alpha * s   */
120:     VecAXPY(z,-alpha,S);   /*     z <- z - alpha * S   */

122:     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
123:       VecNormBegin(r,NORM_2,&dp);
124:     } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
125:       VecNormBegin(z,NORM_2,&dp);
126:     }
127:     VecDotBegin(r,z,&gammaNew);
128:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)r));

130:     KSP_MatMult(ksp,Amat,z,Z);      /*   Z <- Az       */

132:     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
133:       VecNormEnd(r,NORM_2,&dp);
134:     } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
135:       VecNormEnd(z,NORM_2,&dp);
136:     }
137:     VecDotEnd(r,z,&gammaNew);

139:     if (ksp->normtype == KSP_NORM_NATURAL) {
140:       if (PetscIsInfOrNanScalar(gammaNew)) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
141:       dp = PetscSqrtReal(PetscAbsScalar(gammaNew));                  /*     dp <- r'*z = r'*B*r = e'*A'*B*A*e */
142:     } else if (ksp->normtype == KSP_NORM_NONE) {
143:       dp = 0.0;
144:     }
145:     ksp->rnorm = dp;
146:     KSPLogResidualHistory(ksp,dp);
147:     KSPMonitor(ksp,i,dp);
148:     (*ksp->converged)(ksp,i,dp,&ksp->reason,ksp->cnvP);
149:     if (ksp->reason) break;

151:     beta  = gammaNew / gamma;
152:     gamma = gammaNew;
153:     VecAYPX(p,beta,z);   /*     p <- z + beta * p   */
154:     VecAYPX(s,beta,Z);   /*     s <- Z + beta * s   */

156:   } while (i<ksp->max_it);

158:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
159:   return(0);
160: }

164: PETSC_EXTERN PetscErrorCode KSPCreate_GROPPCG(KSP ksp)
165: {

169:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
170:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,1);
171:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
172:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

174:   ksp->ops->setup          = KSPSetUp_GROPPCG;
175:   ksp->ops->solve          = KSPSolve_GROPPCG;
176:   ksp->ops->destroy        = KSPDestroyDefault;
177:   ksp->ops->view           = 0;
178:   ksp->ops->setfromoptions = 0;
179:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
180:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
181:   return(0);
182: }