Actual source code: groppcg.c

petsc-3.4.4 2014-03-13
  1: #include <petsc-private/kspimpl.h>

  3: /*
  4:  KSPSetUp_GROPPCG - Sets up the workspace needed by the GROPPCG method.

  6:  This is called once, usually automatically by KSPSolve() or KSPSetUp()
  7:  but can be called directly by KSPSetUp()
  8: */
 11: PetscErrorCode KSPSetUp_GROPPCG(KSP ksp)
 12: {

 16:   KSPSetWorkVecs(ksp,6);
 17:   return(0);
 18: }

 20: /*
 21:  KSPSolve_GROPPCG

 23:  Input Parameter:
 24:  .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
 25:              example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
 26: */
 29: PetscErrorCode  KSPSolve_GROPPCG(KSP ksp)
 30: {
 32:   PetscInt       i;
 33:   PetscScalar    alpha,beta = 0.0,gamma,gammaNew,t;
 34:   PetscReal      dp = 0.0;
 35:   Vec            x,b,r,p,s,S,z,Z;
 36:   Mat            Amat,Pmat;
 37:   MatStructure   pflag;
 38:   PetscBool      diagonalscale;

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

 44:   x = ksp->vec_sol;
 45:   b = ksp->vec_rhs;
 46:   r = ksp->work[0];
 47:   p = ksp->work[1];
 48:   s = ksp->work[2];
 49:   S = ksp->work[3];
 50:   z = ksp->work[4];
 51:   Z = ksp->work[5];

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

 55:   ksp->its = 0;
 56:   if (!ksp->guess_zero) {
 57:     KSP_MatMult(ksp,Amat,x,r);            /*     r <- b - Ax     */
 58:     VecAYPX(r,-1.0,b);
 59:   } else {
 60:     VecCopy(b,r);                         /*     r <- b (x is 0) */
 61:   }

 63:   KSP_PCApply(ksp,r,z);                   /*     z <- Br   */
 64:   VecCopy(z,p);                           /*     p <- z    */
 65:   VecDotBegin(r,z,&gamma);                  /*     gamma <- z'*r       */
 66:   PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)r));
 67:   KSP_MatMult(ksp,Amat,p,s);              /*     s <- Ap   */
 68:   VecDotEnd(r,z,&gamma);                  /*     gamma <- z'*r       */

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

 94:   i = 0;
 95:   do {
 96:     ksp->its = i+1;
 97:     i++;

 99:     VecDotBegin(p,s,&t);
100:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)p));

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

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

106:     alpha = gamma / t;
107:     VecAXPY(x, alpha,p);   /*     x <- x + alpha * p   */
108:     VecAXPY(r,-alpha,s);   /*     r <- r - alpha * s   */
109:     VecAXPY(z,-alpha,S);   /*     z <- z - alpha * S   */

111:     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
112:       VecNormBegin(r,NORM_2,&dp);
113:     } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
114:       VecNormBegin(z,NORM_2,&dp);
115:     }
116:     VecDotBegin(r,z,&gammaNew);
117:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)r));

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

121:     if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
122:       VecNormEnd(r,NORM_2,&dp);
123:     } else if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
124:       VecNormEnd(z,NORM_2,&dp);
125:     }
126:     VecDotEnd(r,z,&gammaNew);

128:     if (ksp->normtype == KSP_NORM_NATURAL) {
129:       if (PetscIsInfOrNanScalar(gammaNew)) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
130:       dp = PetscSqrtReal(PetscAbsScalar(gammaNew));                  /*     dp <- r'*z = r'*B*r = e'*A'*B*A*e */
131:     } else if (ksp->normtype == KSP_NORM_NONE) {
132:       dp = 0.0;
133:     }
134:     ksp->rnorm = dp;
135:     KSPLogResidualHistory(ksp,dp);
136:     KSPMonitor(ksp,i,dp);
137:     (*ksp->converged)(ksp,i,dp,&ksp->reason,ksp->cnvP);
138:     if (ksp->reason) break;

140:     beta  = gammaNew / gamma;
141:     gamma = gammaNew;
142:     VecAYPX(p,beta,z);   /*     p <- z + beta * p   */
143:     VecAYPX(s,beta,Z);   /*     s <- Z + beta * s   */

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

147:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
148:   return(0);
149: }

151: /*MC
152:    KSPGROPPCG - A pipelined conjugate gradient method from Bill Gropp

154:    There method has two reductions, one of which is overlapped with the matrix-vector product and one of which is
155:    overlapped with the preconditioner.

157:    See also KSPPIPECG, which has only a single reduction that overlaps both the matrix-vector product and the preconditioner.

159:    Level:
160:    beginner

162:    Notes:
163:    MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods.
164:    See the FAQ on the PETSc website for details.

166:    Contributed by:
167:    Pieter Ghysels, Universiteit Antwerpen, Intel Exascience lab Flanders

169:    Reference:
170:    http://www.cs.uiuc.edu/~wgropp/bib/talks/tdata/2012/icerm.pdf

172: .seealso: KSPCreate(), KSPSetType(), KSPPIPECG, KSPPIPECR, KSPPGMRES, KSPCG, KSPCGUseSingleReduction()
173: M*/

177: PETSC_EXTERN PetscErrorCode KSPCreate_GROPPCG(KSP ksp)
178: {

182:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,1);
183:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,1);
184:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,1);
185:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

187:   ksp->ops->setup          = KSPSetUp_GROPPCG;
188:   ksp->ops->solve          = KSPSolve_GROPPCG;
189:   ksp->ops->destroy        = KSPDestroyDefault;
190:   ksp->ops->view           = 0;
191:   ksp->ops->setfromoptions = 0;
192:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
193:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
194:   return(0);
195: }