Actual source code: pipecg.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

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

  4: /*
  5:      KSPSetUp_PIPECG - Sets up the workspace needed by the PIPECG method.

  7:       This is called once, usually automatically by KSPSolve() or KSPSetUp()
  8:      but can be called directly by KSPSetUp()
  9: */
 10: static PetscErrorCode KSPSetUp_PIPECG(KSP ksp)
 11: {

 15:   /* get work vectors needed by PIPECG */
 16:   KSPSetWorkVecs(ksp,9);
 17:   return(0);
 18: }

 20: /*
 21:  KSPSolve_PIPECG - This routine actually applies the pipelined conjugate gradient method

 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: */
 27: static PetscErrorCode  KSPSolve_PIPECG(KSP ksp)
 28: {
 30:   PetscInt       i;
 31:   PetscScalar    alpha = 0.0,beta = 0.0,gamma = 0.0,gammaold = 0.0,delta = 0.0;
 32:   PetscReal      dp    = 0.0;
 33:   Vec            X,B,Z,P,W,Q,U,M,N,R,S;
 34:   Mat            Amat,Pmat;
 35:   PetscBool      diagonalscale;

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

 41:   X = ksp->vec_sol;
 42:   B = ksp->vec_rhs;
 43:   M = ksp->work[0];
 44:   Z = ksp->work[1];
 45:   P = ksp->work[2];
 46:   N = ksp->work[3];
 47:   W = ksp->work[4];
 48:   Q = ksp->work[5];
 49:   U = ksp->work[6];
 50:   R = ksp->work[7];
 51:   S = ksp->work[8];

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

 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,U);                   /*     u <- Br   */

 65:   switch (ksp->normtype) {
 66:   case KSP_NORM_PRECONDITIONED:
 67:     VecNormBegin(U,NORM_2,&dp);                /*     dp <- u'*u = e'*A'*B'*B*A'*e'     */
 68:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)U));
 69:     KSP_MatMult(ksp,Amat,U,W);              /*     w <- Au   */
 70:     VecNormEnd(U,NORM_2,&dp);
 71:     break;
 72:   case KSP_NORM_UNPRECONDITIONED:
 73:     VecNormBegin(R,NORM_2,&dp);                /*     dp <- r'*r = e'*A'*A*e            */
 74:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));
 75:     KSP_MatMult(ksp,Amat,U,W);              /*     w <- Au   */
 76:     VecNormEnd(R,NORM_2,&dp);
 77:     break;
 78:   case KSP_NORM_NATURAL:
 79:     VecDotBegin(R,U,&gamma);                  /*     gamma <- u'*r       */
 80:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));
 81:     KSP_MatMult(ksp,Amat,U,W);              /*     w <- Au   */
 82:     VecDotEnd(R,U,&gamma);
 83:     KSPCheckDot(ksp,gamma);
 84:     dp = PetscSqrtReal(PetscAbsScalar(gamma));                  /*     dp <- r'*u = r'*B*r = e'*A'*B*A*e */
 85:     break;
 86:   case KSP_NORM_NONE:
 87:     KSP_MatMult(ksp,Amat,U,W);
 88:     dp   = 0.0;
 89:     break;
 90:   default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
 91:   }
 92:   KSPLogResidualHistory(ksp,dp);
 93:   KSPMonitor(ksp,0,dp);
 94:   ksp->rnorm = dp;
 95:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
 96:   if (ksp->reason) return(0);

 98:   i = 0;
 99:   do {
100:     if (i > 0 && ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
101:       VecNormBegin(R,NORM_2,&dp);
102:     } else if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
103:       VecNormBegin(U,NORM_2,&dp);
104:     }
105:     if (!(i == 0 && ksp->normtype == KSP_NORM_NATURAL)) {
106:       VecDotBegin(R,U,&gamma);
107:     }
108:     VecDotBegin(W,U,&delta);
109:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));

111:     KSP_PCApply(ksp,W,M);           /*   m <- Bw       */
112:     KSP_MatMult(ksp,Amat,M,N);      /*   n <- Am       */

114:     if (i > 0 && ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
115:       VecNormEnd(R,NORM_2,&dp);
116:     } else if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
117:       VecNormEnd(U,NORM_2,&dp);
118:     }
119:     if (!(i == 0 && ksp->normtype == KSP_NORM_NATURAL)) {
120:       VecDotEnd(R,U,&gamma);
121:     }
122:     VecDotEnd(W,U,&delta);

124:     if (i > 0) {
125:       if (ksp->normtype == KSP_NORM_NATURAL) dp = PetscSqrtReal(PetscAbsScalar(gamma));
126:       else if (ksp->normtype == KSP_NORM_NONE) dp = 0.0;

128:       ksp->rnorm = dp;
129:       KSPLogResidualHistory(ksp,dp);
130:       KSPMonitor(ksp,i,dp);
131:       (*ksp->converged)(ksp,i,dp,&ksp->reason,ksp->cnvP);
132:       if (ksp->reason) break;
133:     }

135:     if (i == 0) {
136:       alpha = gamma / delta;
137:       VecCopy(N,Z);        /*     z <- n          */
138:       VecCopy(M,Q);        /*     q <- m          */
139:       VecCopy(U,P);        /*     p <- u          */
140:       VecCopy(W,S);        /*     s <- w          */
141:     } else {
142:       beta  = gamma / gammaold;
143:       alpha = gamma / (delta - beta / alpha * gamma);
144:       VecAYPX(Z,beta,N);   /*     z <- n + beta * z   */
145:       VecAYPX(Q,beta,M);   /*     q <- m + beta * q   */
146:       VecAYPX(P,beta,U);   /*     p <- u + beta * p   */
147:       VecAYPX(S,beta,W);   /*     s <- w + beta * s   */
148:     }
149:     VecAXPY(X, alpha,P); /*     x <- x + alpha * p   */
150:     VecAXPY(U,-alpha,Q); /*     u <- u - alpha * q   */
151:     VecAXPY(W,-alpha,Z); /*     w <- w - alpha * z   */
152:     VecAXPY(R,-alpha,S); /*     r <- r - alpha * s   */
153:     gammaold = gamma;
154:     i++;
155:     ksp->its = i;

157:     /* if (i%50 == 0) { */
158:     /*   KSP_MatMult(ksp,Amat,X,R);            /\*     w <- b - Ax     *\/ */
159:     /*   VecAYPX(R,-1.0,B); */
160:     /*   KSP_PCApply(ksp,R,U); */
161:     /*   KSP_MatMult(ksp,Amat,U,W); */
162:     /* } */

164:   } while (i<ksp->max_it);
165:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
166:   return(0);
167: }


170: /*MC
171:    KSPPIPECG - Pipelined conjugate gradient method.

173:    This method has only a single non-blocking reduction per iteration, compared to 2 blocking for standard CG.  The
174:    non-blocking reduction is overlapped by the matrix-vector product and preconditioner application.

176:    See also KSPPIPECR, where the reduction is only overlapped with the matrix-vector product.

178:    Level: intermediate

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

184:    Contributed by:
185:    Pieter Ghysels, Universiteit Antwerpen, Intel Exascience lab Flanders

187:    Reference:
188:    P. Ghysels and W. Vanroose, "Hiding global synchronization latency in the preconditioned Conjugate Gradient algorithm",
189:    Submitted to Parallel Computing, 2012.

191: .seealso: KSPCreate(), KSPSetType(), KSPPIPECR, KSPGROPPCG, KSPPGMRES, KSPCG, KSPCGUseSingleReduction()
192: M*/
193: PETSC_EXTERN PetscErrorCode KSPCreate_PIPECG(KSP ksp)
194: {

198:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
199:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
200:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
201:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

203:   ksp->ops->setup          = KSPSetUp_PIPECG;
204:   ksp->ops->solve          = KSPSolve_PIPECG;
205:   ksp->ops->destroy        = KSPDestroyDefault;
206:   ksp->ops->view           = 0;
207:   ksp->ops->setfromoptions = 0;
208:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
209:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
210:   return(0);
211: }