Actual source code: pipecg.c
petsc-3.7.5 2017-01-01
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: */
12: PetscErrorCode KSPSetUp_PIPECG(KSP ksp)
13: {
17: /* get work vectors needed by PIPECG */
18: KSPSetWorkVecs(ksp,9);
19: return(0);
20: }
22: /*
23: KSPSolve_PIPECG - This routine actually applies the pipelined conjugate gradient method
25: Input Parameter:
26: . ksp - the Krylov space object that was set to use conjugate gradient, by, for
27: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
28: */
31: PetscErrorCode KSPSolve_PIPECG(KSP ksp)
32: {
34: PetscInt i;
35: PetscScalar alpha = 0.0,beta = 0.0,gamma = 0.0,gammaold = 0.0,delta = 0.0;
36: PetscReal dp = 0.0;
37: Vec X,B,Z,P,W,Q,U,M,N,R,S;
38: Mat Amat,Pmat;
39: PetscBool diagonalscale;
42: PCGetDiagonalScale(ksp->pc,&diagonalscale);
43: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
45: X = ksp->vec_sol;
46: B = ksp->vec_rhs;
47: M = ksp->work[0];
48: Z = ksp->work[1];
49: P = ksp->work[2];
50: N = ksp->work[3];
51: W = ksp->work[4];
52: Q = ksp->work[5];
53: U = ksp->work[6];
54: R = ksp->work[7];
55: S = ksp->work[8];
57: PCGetOperators(ksp->pc,&Amat,&Pmat);
59: ksp->its = 0;
60: if (!ksp->guess_zero) {
61: KSP_MatMult(ksp,Amat,X,R); /* r <- b - Ax */
62: VecAYPX(R,-1.0,B);
63: } else {
64: VecCopy(B,R); /* r <- b (x is 0) */
65: }
67: KSP_PCApply(ksp,R,U); /* u <- Br */
69: switch (ksp->normtype) {
70: case KSP_NORM_PRECONDITIONED:
71: VecNormBegin(U,NORM_2,&dp); /* dp <- u'*u = e'*A'*B'*B*A'*e' */
72: PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)U));
73: KSP_MatMult(ksp,Amat,U,W); /* w <- Au */
74: VecNormEnd(U,NORM_2,&dp);
75: break;
76: case KSP_NORM_UNPRECONDITIONED:
77: VecNormBegin(R,NORM_2,&dp); /* dp <- r'*r = e'*A'*A*e */
78: PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));
79: KSP_MatMult(ksp,Amat,U,W); /* w <- Au */
80: VecNormEnd(R,NORM_2,&dp);
81: break;
82: case KSP_NORM_NATURAL:
83: VecDotBegin(R,U,&gamma); /* gamma <- u'*r */
84: PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));
85: KSP_MatMult(ksp,Amat,U,W); /* w <- Au */
86: VecDotEnd(R,U,&gamma);
87: KSPCheckDot(ksp,gamma);
88: dp = PetscSqrtReal(PetscAbsScalar(gamma)); /* dp <- r'*u = r'*B*r = e'*A'*B*A*e */
89: break;
90: case KSP_NORM_NONE:
91: KSP_MatMult(ksp,Amat,U,W);
92: dp = 0.0;
93: break;
94: default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
95: }
96: KSPLogResidualHistory(ksp,dp);
97: KSPMonitor(ksp,0,dp);
98: ksp->rnorm = dp;
99: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
100: if (ksp->reason) return(0);
102: i = 0;
103: do {
104: if (i > 0 && ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
105: VecNormBegin(R,NORM_2,&dp);
106: } else if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
107: VecNormBegin(U,NORM_2,&dp);
108: }
109: if (!(i == 0 && ksp->normtype == KSP_NORM_NATURAL)) {
110: VecDotBegin(R,U,&gamma);
111: }
112: VecDotBegin(W,U,&delta);
113: PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)R));
115: KSP_PCApply(ksp,W,M); /* m <- Bw */
116: KSP_MatMult(ksp,Amat,M,N); /* n <- Am */
118: if (i > 0 && ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
119: VecNormEnd(R,NORM_2,&dp);
120: } else if (i > 0 && ksp->normtype == KSP_NORM_PRECONDITIONED) {
121: VecNormEnd(U,NORM_2,&dp);
122: }
123: if (!(i == 0 && ksp->normtype == KSP_NORM_NATURAL)) {
124: VecDotEnd(R,U,&gamma);
125: }
126: VecDotEnd(W,U,&delta);
128: if (i > 0) {
129: if (ksp->normtype == KSP_NORM_NATURAL) dp = PetscSqrtReal(PetscAbsScalar(gamma));
130: else if (ksp->normtype == KSP_NORM_NONE) dp = 0.0;
132: ksp->rnorm = dp;
133: KSPLogResidualHistory(ksp,dp);
134: KSPMonitor(ksp,i,dp);
135: (*ksp->converged)(ksp,i,dp,&ksp->reason,ksp->cnvP);
136: if (ksp->reason) break;
137: }
139: if (i == 0) {
140: alpha = gamma / delta;
141: VecCopy(N,Z); /* z <- n */
142: VecCopy(M,Q); /* q <- m */
143: VecCopy(U,P); /* p <- u */
144: VecCopy(W,S); /* s <- w */
145: } else {
146: beta = gamma / gammaold;
147: alpha = gamma / (delta - beta / alpha * gamma);
148: VecAYPX(Z,beta,N); /* z <- n + beta * z */
149: VecAYPX(Q,beta,M); /* q <- m + beta * q */
150: VecAYPX(P,beta,U); /* p <- u + beta * p */
151: VecAYPX(S,beta,W); /* s <- w + beta * s */
152: }
153: VecAXPY(X, alpha,P); /* x <- x + alpha * p */
154: VecAXPY(U,-alpha,Q); /* u <- u - alpha * q */
155: VecAXPY(W,-alpha,Z); /* w <- w - alpha * z */
156: VecAXPY(R,-alpha,S); /* r <- r - alpha * s */
157: gammaold = gamma;
158: i++;
159: ksp->its = i;
161: /* if (i%50 == 0) { */
162: /* KSP_MatMult(ksp,Amat,X,R); /\* w <- b - Ax *\/ */
163: /* VecAYPX(R,-1.0,B); */
164: /* KSP_PCApply(ksp,R,U); */
165: /* KSP_MatMult(ksp,Amat,U,W); */
166: /* } */
168: } while (i<ksp->max_it);
169: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
170: return(0);
171: }
174: /*MC
175: KSPPIPECG - Pipelined conjugate gradient method.
177: This method has only a single non-blocking reduction per iteration, compared to 2 blocking for standard CG. The
178: non-blocking reduction is overlapped by the matrix-vector product and preconditioner application.
180: See also KSPPIPECR, where the reduction is only overlapped with the matrix-vector product.
182: Level: intermediate
184: Notes:
185: MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods.
186: See the FAQ on the PETSc website for details.
188: Contributed by:
189: Pieter Ghysels, Universiteit Antwerpen, Intel Exascience lab Flanders
191: Reference:
192: P. Ghysels and W. Vanroose, "Hiding global synchronization latency in the preconditioned Conjugate Gradient algorithm",
193: Submitted to Parallel Computing, 2012.
195: .seealso: KSPCreate(), KSPSetType(), KSPPIPECR, KSPGROPPCG, KSPPGMRES, KSPCG, KSPCGUseSingleReduction()
196: M*/
199: PETSC_EXTERN PetscErrorCode KSPCreate_PIPECG(KSP ksp)
200: {
204: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
205: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
206: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
208: ksp->ops->setup = KSPSetUp_PIPECG;
209: ksp->ops->solve = KSPSolve_PIPECG;
210: ksp->ops->destroy = KSPDestroyDefault;
211: ksp->ops->view = 0;
212: ksp->ops->setfromoptions = 0;
213: ksp->ops->buildsolution = KSPBuildSolutionDefault;
214: ksp->ops->buildresidual = KSPBuildResidualDefault;
215: return(0);
216: }