Actual source code: pipecr.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: #include <petsc/private/kspimpl.h>

  3: /*
  4:      KSPSetUp_PIPECR - Sets up the workspace needed by the PIPECR method.

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

 16:   /* get work vectors needed by PIPECR */
 17:   KSPSetWorkVecs(ksp,7);
 18:   return(0);
 19: }

 21: /*
 22:  KSPSolve_PIPECR - This routine actually applies the pipelined conjugate residual method

 24:  Input Parameter:
 25:  .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
 26:              example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
 27: */
 30: PetscErrorCode  KSPSolve_PIPECR(KSP ksp)
 31: {
 33:   PetscInt       i;
 34:   PetscScalar    alpha=0.0,beta=0.0,gamma,gammaold=0.0,delta;
 35:   PetscReal      dp   = 0.0;
 36:   Vec            X,B,Z,P,W,Q,U,M,N;
 37:   Mat            Amat,Pmat;
 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:   M = ksp->work[0];
 47:   Z = ksp->work[1];
 48:   P = ksp->work[2];
 49:   N = ksp->work[3];
 50:   W = ksp->work[4];
 51:   Q = ksp->work[5];
 52:   U = ksp->work[6];

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

 56:   ksp->its = 0;
 57:   /* we don't have an R vector, so put the (unpreconditioned) residual in w for now */
 58:   if (!ksp->guess_zero) {
 59:     KSP_MatMult(ksp,Amat,X,W);            /*     w <- b - Ax     */
 60:     VecAYPX(W,-1.0,B);
 61:   } else {
 62:     VecCopy(B,W);                         /*     w <- b (x is 0) */
 63:   }
 64:   KSP_PCApply(ksp,W,U);                   /*     u <- Bw   */

 66:   switch (ksp->normtype) {
 67:   case KSP_NORM_PRECONDITIONED:
 68:     VecNormBegin(U,NORM_2,&dp);           /*     dp <- u'*u = e'*A'*B'*B*A'*e'     */
 69:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)U));
 70:     KSP_MatMult(ksp,Amat,U,W);            /*     w <- Au   */
 71:     VecNormEnd(U,NORM_2,&dp);
 72:     break;
 73:   case KSP_NORM_NONE:
 74:     KSP_MatMult(ksp,Amat,U,W);
 75:     dp   = 0.0;
 76:     break;
 77:   default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
 78:   }
 79:   KSPLogResidualHistory(ksp,dp);
 80:   KSPMonitor(ksp,0,dp);
 81:   ksp->rnorm = dp;
 82:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
 83:   if (ksp->reason) return(0);

 85:   i = 0;
 86:   do {
 87:     KSP_PCApply(ksp,W,M);            /*   m <- Bw       */

 89:     if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
 90:       VecNormBegin(U,NORM_2,&dp);
 91:     }
 92:     VecDotBegin(W,U,&gamma);
 93:     VecDotBegin(M,W,&delta);
 94:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)U));

 96:     KSP_MatMult(ksp,Amat,M,N);       /*   n <- Am       */

 98:     if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
 99:       VecNormEnd(U,NORM_2,&dp);
100:     }
101:     VecDotEnd(W,U,&gamma);
102:     VecDotEnd(M,W,&delta);

104:     if (i > 0) {
105:       if (ksp->normtype == KSP_NORM_NONE) dp = 0.0;
106:       ksp->rnorm = dp;
107:       KSPLogResidualHistory(ksp,dp);
108:       KSPMonitor(ksp,i,dp);
109:       (*ksp->converged)(ksp,i,dp,&ksp->reason,ksp->cnvP);
110:       if (ksp->reason) break;
111:     }

113:     if (i == 0) {
114:       alpha = gamma / delta;
115:       VecCopy(N,Z);        /*     z <- n          */
116:       VecCopy(M,Q);        /*     q <- m          */
117:       VecCopy(U,P);        /*     p <- u          */
118:     } else {
119:       beta  = gamma / gammaold;
120:       alpha = gamma / (delta - beta / alpha * gamma);
121:       VecAYPX(Z,beta,N);   /*     z <- n + beta * z   */
122:       VecAYPX(Q,beta,M);   /*     q <- m + beta * q   */
123:       VecAYPX(P,beta,U);   /*     p <- u + beta * p   */
124:     }
125:     VecAXPY(X, alpha,P); /*     x <- x + alpha * p   */
126:     VecAXPY(U,-alpha,Q); /*     u <- u - alpha * q   */
127:     VecAXPY(W,-alpha,Z); /*     w <- w - alpha * z   */
128:     gammaold = gamma;
129:     i++;
130:     ksp->its = i;

132:     /* if (i%50 == 0) { */
133:     /*   KSP_MatMult(ksp,Amat,X,W);            /\*     w <- b - Ax     *\/ */
134:     /*   VecAYPX(W,-1.0,B); */
135:     /*   KSP_PCApply(ksp,W,U); */
136:     /*   KSP_MatMult(ksp,Amat,U,W); */
137:     /* } */

139:   } while (i<ksp->max_it);
140:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
141:   return(0);
142: }

144: /*MC
145:    KSPPIPECR - Pipelined conjugate residual method

147:    This method has only a single non-blocking reduction per iteration, compared to 2 blocking for standard CR.  The
148:    non-blocking reduction is overlapped by the matrix-vector product, but not the preconditioner application.

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

152:    Level: intermediate

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

158:    Contributed by:
159:    Pieter Ghysels, Universiteit Antwerpen, Intel Exascience lab Flanders

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

165: .seealso: KSPCreate(), KSPSetType(), KSPPIPECG, KSPGROPPCG, KSPPGMRES, KSPCG, KSPCGUseSingleReduction()
166: M*/

170: PETSC_EXTERN PetscErrorCode KSPCreate_PIPECR(KSP ksp)
171: {

175:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);

177:   ksp->ops->setup          = KSPSetUp_PIPECR;
178:   ksp->ops->solve          = KSPSolve_PIPECR;
179:   ksp->ops->destroy        = KSPDestroyDefault;
180:   ksp->ops->view           = 0;
181:   ksp->ops->setfromoptions = 0;
182:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
183:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
184:   return(0);
185: }