Actual source code: pipecg.c


  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: {
 12:   /* get work vectors needed by PIPECG */
 13:   KSPSetWorkVecs(ksp,9);
 14:   return 0;
 15: }

 17: /*
 18:  KSPSolve_PIPECG - This routine actually applies the pipelined conjugate gradient method
 19: */
 20: static PetscErrorCode  KSPSolve_PIPECG(KSP ksp)
 21: {
 22:   PetscInt       i;
 23:   PetscScalar    alpha = 0.0,beta = 0.0,gamma = 0.0,gammaold = 0.0,delta = 0.0;
 24:   PetscReal      dp    = 0.0;
 25:   Vec            X,B,Z,P,W,Q,U,M,N,R,S;
 26:   Mat            Amat,Pmat;
 27:   PetscBool      diagonalscale;

 29:   PCGetDiagonalScale(ksp->pc,&diagonalscale);

 32:   X = ksp->vec_sol;
 33:   B = ksp->vec_rhs;
 34:   R = ksp->work[0];
 35:   Z = ksp->work[1];
 36:   P = ksp->work[2];
 37:   N = ksp->work[3];
 38:   W = ksp->work[4];
 39:   Q = ksp->work[5];
 40:   U = ksp->work[6];
 41:   M = ksp->work[7];
 42:   S = ksp->work[8];

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

 46:   ksp->its = 0;
 47:   if (!ksp->guess_zero) {
 48:     KSP_MatMult(ksp,Amat,X,R);            /*     r <- b - Ax     */
 49:     VecAYPX(R,-1.0,B);
 50:   } else {
 51:     VecCopy(B,R);                         /*     r <- b (x is 0) */
 52:   }

 54:   KSP_PCApply(ksp,R,U);                   /*     u <- Br   */

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

 89:   i = 0;
 90:   do {
 91:     if (i > 0 && ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
 92:       VecNormBegin(R,NORM_2,&dp);
 93:     } else if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
 94:       VecNormBegin(U,NORM_2,&dp);
 95:     }
 96:     if (!(i == 0 && ksp->normtype == KSP_NORM_NATURAL)) {
 97:       VecDotBegin(R,U,&gamma);
 98:     }
 99:     VecDotBegin(W,U,&delta);
100:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));

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

105:     if (i > 0 && ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
106:       VecNormEnd(R,NORM_2,&dp);
107:     } else if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
108:       VecNormEnd(U,NORM_2,&dp);
109:     }
110:     if (!(i == 0 && ksp->normtype == KSP_NORM_NATURAL)) {
111:       VecDotEnd(R,U,&gamma);
112:     }
113:     VecDotEnd(W,U,&delta);

115:     if (i > 0) {
116:       if (ksp->normtype == KSP_NORM_NATURAL) dp = PetscSqrtReal(PetscAbsScalar(gamma));
117:       else if (ksp->normtype == KSP_NORM_NONE) dp = 0.0;

119:       ksp->rnorm = dp;
120:       KSPLogResidualHistory(ksp,dp);
121:       KSPMonitor(ksp,i,dp);
122:       (*ksp->converged)(ksp,i,dp,&ksp->reason,ksp->cnvP);
123:       if (ksp->reason) return 0;
124:     }

126:     if (i == 0) {
127:       alpha = gamma / delta;
128:       VecCopy(N,Z);        /*     z <- n          */
129:       VecCopy(M,Q);        /*     q <- m          */
130:       VecCopy(U,P);        /*     p <- u          */
131:       VecCopy(W,S);        /*     s <- w          */
132:     } else {
133:       beta  = gamma / gammaold;
134:       alpha = gamma / (delta - beta / alpha * gamma);
135:       VecAYPX(Z,beta,N);   /*     z <- n + beta * z   */
136:       VecAYPX(Q,beta,M);   /*     q <- m + beta * q   */
137:       VecAYPX(P,beta,U);   /*     p <- u + beta * p   */
138:       VecAYPX(S,beta,W);   /*     s <- w + beta * s   */
139:     }
140:     VecAXPY(X, alpha,P); /*     x <- x + alpha * p   */
141:     VecAXPY(U,-alpha,Q); /*     u <- u - alpha * q   */
142:     VecAXPY(W,-alpha,Z); /*     w <- w - alpha * z   */
143:     VecAXPY(R,-alpha,S); /*     r <- r - alpha * s   */
144:     gammaold = gamma;
145:     i++;
146:     ksp->its = i;

148:     /* if (i%50 == 0) { */
149:     /*   KSP_MatMult(ksp,Amat,X,R);            /\*     w <- b - Ax     *\/ */
150:     /*   VecAYPX(R,-1.0,B); */
151:     /*   KSP_PCApply(ksp,R,U); */
152:     /*   KSP_MatMult(ksp,Amat,U,W); */
153:     /* } */

155:   } while (i<=ksp->max_it);
156:   if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
157:   return 0;
158: }

160: PETSC_INTERN PetscErrorCode KSPBuildResidual_CG(KSP,Vec,Vec,Vec*);

162: /*MC
163:    KSPPIPECG - Pipelined conjugate gradient method.

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

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

170:    Level: intermediate

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

176:    Contributed by:
177:    Pieter Ghysels, Universiteit Antwerpen, Intel Exascience lab Flanders

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

183: .seealso: KSPCreate(), KSPSetType(), KSPPIPECR, KSPGROPPCG, KSPPGMRES, KSPCG, KSPCGUseSingleReduction()
184: M*/
185: PETSC_EXTERN PetscErrorCode KSPCreate_PIPECG(KSP ksp)
186: {
187:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
188:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
189:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
190:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

192:   ksp->ops->setup          = KSPSetUp_PIPECG;
193:   ksp->ops->solve          = KSPSolve_PIPECG;
194:   ksp->ops->destroy        = KSPDestroyDefault;
195:   ksp->ops->view           = NULL;
196:   ksp->ops->setfromoptions = NULL;
197:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
198:   ksp->ops->buildresidual  = KSPBuildResidual_CG;
199:   return 0;
200: }