xref: /petsc/src/ksp/ksp/impls/cg/cg.c (revision 65f8aed5f7eaa1e2ef2ddeffe666264e0669c876)
1 
2 /*
3     This file implements the conjugate gradient method in PETSc as part of
4     KSP. You can use this as a starting point for implementing your own
5     Krylov method that is not provided with PETSc.
6 
7     The following basic routines are required for each Krylov method.
8         KSPCreate_XXX()          - Creates the Krylov context
9         KSPSetFromOptions_XXX()  - Sets runtime options
10         KSPSolve_XXX()           - Runs the Krylov method
11         KSPDestroy_XXX()         - Destroys the Krylov context, freeing all
12                                    memory it needed
13     Here the "_XXX" denotes a particular implementation, in this case
14     we use _CG (e.g. KSPCreate_CG, KSPDestroy_CG). These routines are
15     are actually called via the common user interface routines
16     KSPSetType(), KSPSetFromOptions(), KSPSolve(), and KSPDestroy() so the
17     application code interface remains identical for all preconditioners.
18 
19     Other basic routines for the KSP objects include
20         KSPSetUp_XXX()
21         KSPView_XXX()            - Prints details of solver being used.
22 
23     Detailed Notes:
24     By default, this code implements the CG (Conjugate Gradient) method,
25     which is valid for real symmetric (and complex Hermitian) positive
26     definite matrices. Note that for the complex Hermitian case, the
27     VecDot() arguments within the code MUST remain in the order given
28     for correct computation of inner products.
29 
30     Reference: Hestenes and Steifel, 1952.
31 
32     By switching to the indefinite vector inner product, VecTDot(), the
33     same code is used for the complex symmetric case as well.  The user
34     must call KSPCGSetType(ksp,KSP_CG_SYMMETRIC) or use the option
35     -ksp_cg_type symmetric to invoke this variant for the complex case.
36     Note, however, that the complex symmetric code is NOT valid for
37     all such matrices ... and thus we don't recommend using this method.
38 */
39 /*
40     cgimpl.h defines the simple data structured used to store information
41     related to the type of matrix (e.g. complex symmetric) being solved and
42     data used during the optional Lanczo process used to compute eigenvalues
43 */
44 #include <../src/ksp/ksp/impls/cg/cgimpl.h>       /*I "petscksp.h" I*/
45 extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal*,PetscReal*);
46 extern PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*);
47 
48 /*
49      KSPSetUp_CG - Sets up the workspace needed by the CG method.
50 
51       This is called once, usually automatically by KSPSolve() or KSPSetUp()
52      but can be called directly by KSPSetUp()
53 */
54 static PetscErrorCode KSPSetUp_CG(KSP ksp)
55 {
56   KSP_CG         *cgP = (KSP_CG*)ksp->data;
57   PetscErrorCode ierr;
58   PetscInt       maxit = ksp->max_it,nwork = 3;
59 
60   PetscFunctionBegin;
61   /* get work vectors needed by CG */
62   if (cgP->singlereduction) nwork += 2;
63   ierr = KSPSetWorkVecs(ksp,nwork);CHKERRQ(ierr);
64 
65   /*
66      If user requested computations of eigenvalues then allocate work
67      work space needed
68   */
69   if (ksp->calc_sings) {
70     /* get space to store tridiagonal matrix for Lanczos */
71     ierr = PetscMalloc4(maxit+1,&cgP->e,maxit+1,&cgP->d,maxit+1,&cgP->ee,maxit+1,&cgP->dd);CHKERRQ(ierr);
72     ierr = PetscLogObjectMemory((PetscObject)ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));CHKERRQ(ierr);
73 
74     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
75     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
76   }
77   PetscFunctionReturn(0);
78 }
79 
80 /*
81      A macro used in the following KSPSolve_CG and KSPSolve_CG_SingleReduction routines
82 */
83 #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
84 
85 /*
86      KSPSolve_CG - This routine actually applies the conjugate gradient method
87 
88      Note : this routine can be replaced with another one (see below) which implements
89             another variant of CG.
90 
91    Input Parameter:
92 .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
93             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
94 */
95 static PetscErrorCode KSPSolve_CG(KSP ksp)
96 {
97   PetscErrorCode ierr;
98   PetscInt       i,stored_max_it,eigs;
99   PetscScalar    dpi = 0.0,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0,dpiold;
100   PetscReal      dp  = 0.0;
101   Vec            X,B,Z,R,P,W;
102   KSP_CG         *cg;
103   Mat            Amat,Pmat;
104   PetscBool      diagonalscale;
105 
106   PetscFunctionBegin;
107   ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr);
108   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
109 
110   cg            = (KSP_CG*)ksp->data;
111   eigs          = ksp->calc_sings;
112   stored_max_it = ksp->max_it;
113   X             = ksp->vec_sol;
114   B             = ksp->vec_rhs;
115   R             = ksp->work[0];
116   Z             = ksp->work[1];
117   P             = ksp->work[2];
118   W             = Z;
119 
120   if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
121   ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr);
122 
123   ksp->its = 0;
124   if (!ksp->guess_zero) {
125     ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr);            /*    r <- b - Ax                       */
126     ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr);
127   } else {
128     ierr = VecCopy(B,R);CHKERRQ(ierr);                         /*    r <- b (x is 0)                   */
129   }
130 
131   switch (ksp->normtype) {
132     case KSP_NORM_PRECONDITIONED:
133       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
134       ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr);              /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
135       break;
136     case KSP_NORM_UNPRECONDITIONED:
137       ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr);              /*    dp <- r'*r = e'*A'*A*e            */
138       break;
139     case KSP_NORM_NATURAL:
140       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
141       ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr);                 /*    beta <- z'*r                      */
142       KSPCheckDot(ksp,beta);
143       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
144       break;
145     case KSP_NORM_NONE:
146       dp = 0.0;
147       break;
148     default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
149   }
150   ierr       = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr);
151   ierr       = KSPMonitor(ksp,0,dp);CHKERRQ(ierr);
152   ksp->rnorm = dp;
153 
154   ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);     /* test for convergence */
155   if (ksp->reason) PetscFunctionReturn(0);
156 
157   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) {
158     ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);                /*     z <- Br                           */
159   }
160   if (ksp->normtype != KSP_NORM_NATURAL) {
161     ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr);                  /*     beta <- z'*r                      */
162     KSPCheckDot(ksp,beta);
163   }
164 
165   i = 0;
166   do {
167     ksp->its = i+1;
168     if (beta == 0.0) {
169       ksp->reason = KSP_CONVERGED_ATOL;
170       ierr        = PetscInfo(ksp,"converged due to beta = 0\n");CHKERRQ(ierr);
171       break;
172 #if !defined(PETSC_USE_COMPLEX)
173     } else if ((i > 0) && (beta*betaold < 0.0)) {
174       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite preconditioner");
175       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
176       ierr        = PetscInfo(ksp,"diverging due to indefinite preconditioner\n");CHKERRQ(ierr);
177       break;
178 #endif
179     }
180     if (!i) {
181       ierr = VecCopy(Z,P);CHKERRQ(ierr);                       /*     p <- z                           */
182       b    = 0.0;
183     } else {
184       b = beta/betaold;
185       if (eigs) {
186         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
187         e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
188       }
189       ierr = VecAYPX(P,b,Z);CHKERRQ(ierr);                     /*     p <- z + b* p                    */
190     }
191     dpiold = dpi;
192     ierr = KSP_MatMult(ksp,Amat,P,W);CHKERRQ(ierr);            /*     w <- Ap                          */
193     ierr = VecXDot(P,W,&dpi);CHKERRQ(ierr);                    /*     dpi <- p'w                       */
194     KSPCheckDot(ksp,dpi);
195     betaold = beta;
196 
197     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi*dpiold) <= 0.0))) {
198       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite matrix");
199       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
200       ierr        = PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");CHKERRQ(ierr);
201       break;
202     }
203     a = beta/dpi;                                              /*     a = beta/p'w                     */
204     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
205     ierr = VecAXPY(X,a,P);CHKERRQ(ierr);                       /*     x <- x + ap                      */
206     ierr = VecAXPY(R,-a,W);CHKERRQ(ierr);                      /*     r <- r - aw                      */
207     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
208       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*     z <- Br                          */
209       ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr);              /*     dp <- z'*z                       */
210     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
211       ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr);              /*     dp <- r'*r                       */
212     } else if (ksp->normtype == KSP_NORM_NATURAL) {
213       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*     z <- Br                          */
214       ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr);                 /*     beta <- r'*z                     */
215       KSPCheckDot(ksp,beta);
216       dp = PetscSqrtReal(PetscAbsScalar(beta));
217     } else {
218       dp = 0.0;
219     }
220     ksp->rnorm = dp;
221     ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr);
222     if (eigs) cg->ned = ksp->its;
223     ierr = KSPMonitor(ksp,i+1,dp);CHKERRQ(ierr);
224     ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
225     if (ksp->reason) break;
226 
227     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) {
228       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*     z <- Br                          */
229     }
230     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) {
231       ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr);                 /*     beta <- z'*r                     */
232       KSPCheckDot(ksp,beta);
233     }
234 
235     i++;
236   } while (i<ksp->max_it);
237   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
238   PetscFunctionReturn(0);
239 }
240 
241 /*
242        KSPSolve_CG_SingleReduction
243 
244        This variant of CG is identical in exact arithmetic to the standard algorithm,
245        but is rearranged to use only a single reduction stage per iteration, using additional
246        intermediate vectors.
247 
248        See KSPCGUseSingleReduction_CG()
249 
250 */
251 static PetscErrorCode KSPSolve_CG_SingleReduction(KSP ksp)
252 {
253   PetscErrorCode ierr;
254   PetscInt       i,stored_max_it,eigs;
255   PetscScalar    dpi = 0.0,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0,delta,dpiold,tmp[2];
256   PetscReal      dp  = 0.0;
257   Vec            X,B,Z,R,P,S,W,tmpvecs[2];
258   KSP_CG         *cg;
259   Mat            Amat,Pmat;
260   PetscBool      diagonalscale;
261 
262   PetscFunctionBegin;
263   ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr);
264   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
265 
266   cg            = (KSP_CG*)ksp->data;
267   eigs          = ksp->calc_sings;
268   stored_max_it = ksp->max_it;
269   X             = ksp->vec_sol;
270   B             = ksp->vec_rhs;
271   R             = ksp->work[0];
272   Z             = ksp->work[1];
273   P             = ksp->work[2];
274   S             = ksp->work[3];
275   W             = ksp->work[4];
276 
277   if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
278   ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr);
279 
280   ksp->its = 0;
281   if (!ksp->guess_zero) {
282     ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr);            /*    r <- b - Ax                       */
283     ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr);
284   } else {
285     ierr = VecCopy(B,R);CHKERRQ(ierr);                         /*    r <- b (x is 0)                   */
286   }
287 
288   switch (ksp->normtype) {
289     case KSP_NORM_PRECONDITIONED:
290       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
291       ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr);              /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
292       break;
293     case KSP_NORM_UNPRECONDITIONED:
294       ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr);              /*    dp <- r'*r = e'*A'*A*e            */
295       break;
296     case KSP_NORM_NATURAL:
297       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
298       ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr);
299       ierr = VecXDot(Z,S,&delta);CHKERRQ(ierr);                /*    delta <- z'*A*z = r'*B*A*B*r      */
300       ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr);                 /*    beta <- z'*r                      */
301       KSPCheckDot(ksp,beta);
302       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
303       break;
304     case KSP_NORM_NONE:
305       dp = 0.0;
306       break;
307     default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]);
308   }
309   ierr       = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr);
310   ierr       = KSPMonitor(ksp,0,dp);CHKERRQ(ierr);
311   ksp->rnorm = dp;
312 
313   ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);     /* test for convergence */
314   if (ksp->reason) PetscFunctionReturn(0);
315 
316   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) {
317     ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);                 /*    z <- Br                           */
318   }
319   if (ksp->normtype != KSP_NORM_NATURAL) {
320     ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr);
321     ierr = VecXDot(Z,S,&delta);CHKERRQ(ierr);                  /*    delta <- z'*A*z = r'*B*A*B*r      */
322     ierr = VecXDot(Z,R,&beta);CHKERRQ(ierr);                   /*    beta <- z'*r                      */
323     KSPCheckDot(ksp,beta);
324   }
325 
326   i = 0;
327   do {
328     ksp->its = i+1;
329     if (beta == 0.0) {
330       ksp->reason = KSP_CONVERGED_ATOL;
331       ierr        = PetscInfo(ksp,"converged due to beta = 0\n");CHKERRQ(ierr);
332       break;
333 #if !defined(PETSC_USE_COMPLEX)
334     } else if ((i > 0) && (beta*betaold < 0.0)) {
335       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite preconditioner");
336       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
337       ierr        = PetscInfo(ksp,"diverging due to indefinite preconditioner\n");CHKERRQ(ierr);
338       break;
339 #endif
340     }
341     if (!i) {
342       ierr = VecCopy(Z,P);CHKERRQ(ierr);                       /*    p <- z                           */
343       b    = 0.0;
344     } else {
345       b = beta/betaold;
346       if (eigs) {
347         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
348         e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
349       }
350       ierr = VecAYPX(P,b,Z);CHKERRQ(ierr);                     /*    p <- z + b* p                     */
351     }
352     dpiold = dpi;
353     if (!i) {
354       ierr = KSP_MatMult(ksp,Amat,P,W);CHKERRQ(ierr);          /*    w <- Ap                           */
355       ierr = VecXDot(P,W,&dpi);CHKERRQ(ierr);                  /*    dpi <- p'w                        */
356     } else {
357       ierr = VecAYPX(W,beta/betaold,S);CHKERRQ(ierr);          /*    w <- Ap                           */
358       dpi  = delta - beta*beta*dpiold/(betaold*betaold);       /*    dpi <- p'w                        */
359     }
360     betaold = beta;
361     KSPCheckDot(ksp,beta);
362 
363     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi*dpiold) <= 0.0))) {
364       if (ksp->errorifnotconverged) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"Diverged due to indefinite matrix");
365       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
366       ierr        = PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");CHKERRQ(ierr);
367       break;
368     }
369     a = beta/dpi;                                              /*    a = beta/p'w                      */
370     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
371     ierr = VecAXPY(X,a,P);CHKERRQ(ierr);                       /*    x <- x + ap                       */
372     ierr = VecAXPY(R,-a,W);CHKERRQ(ierr);                      /*    r <- r - aw                       */
373     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
374       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
375       ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr);
376       ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr);              /*    dp <- z'*z                        */
377     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
378       ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr);              /*    dp <- r'*r                        */
379     } else if (ksp->normtype == KSP_NORM_NATURAL) {
380       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
381       tmpvecs[0] = S; tmpvecs[1] = R;
382       ierr  = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr);
383       ierr  = VecMDot(Z,2,tmpvecs,tmp);CHKERRQ(ierr);          /*    delta <- z'*A*z = r'*B*A*B*r      */
384       delta = tmp[0]; beta = tmp[1];                           /*    beta <- z'*r                      */
385       KSPCheckDot(ksp,beta);
386       dp = PetscSqrtReal(PetscAbsScalar(beta));                /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
387     } else {
388       dp = 0.0;
389     }
390     ksp->rnorm = dp;
391     CHKERRQ(ierr);KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr);
392     if (eigs) cg->ned = ksp->its;
393     ierr = KSPMonitor(ksp,i+1,dp);CHKERRQ(ierr);
394     ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
395     if (ksp->reason) break;
396 
397     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)) {
398       ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr);               /*    z <- Br                           */
399       ierr = KSP_MatMult(ksp,Amat,Z,S);CHKERRQ(ierr);
400     }
401     if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)) {
402       tmpvecs[0] = S; tmpvecs[1] = R;
403       ierr  = VecMDot(Z,2,tmpvecs,tmp);CHKERRQ(ierr);
404       delta = tmp[0]; beta = tmp[1];                           /*    delta <- z'*A*z = r'*B'*A*B*r     */
405       KSPCheckDot(ksp,beta);                                   /*    beta <- z'*r                      */
406     }
407 
408     i++;
409   } while (i<ksp->max_it);
410   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
411   PetscFunctionReturn(0);
412 }
413 
414 /*
415      KSPDestroy_CG - Frees resources allocated in KSPSetup_CG and clears function
416                      compositions from KSPCreate_CG. If adding your own KSP implementation,
417                      you must be sure to free all allocated resources here to prevent
418                      leaks.
419 */
420 PetscErrorCode KSPDestroy_CG(KSP ksp)
421 {
422   KSP_CG         *cg = (KSP_CG*)ksp->data;
423   PetscErrorCode ierr;
424 
425   PetscFunctionBegin;
426   /* free space used for singular value calculations */
427   if (ksp->calc_sings) {
428     ierr = PetscFree4(cg->e,cg->d,cg->ee,cg->dd);CHKERRQ(ierr);
429   }
430   ierr = KSPDestroyDefault(ksp);CHKERRQ(ierr);
431   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",NULL);CHKERRQ(ierr);
432   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",NULL);CHKERRQ(ierr);
433   PetscFunctionReturn(0);
434 }
435 
436 /*
437      KSPView_CG - Prints information about the current Krylov method being used.
438                   If your Krylov method has special options or flags that information
439                   should be printed here.
440 */
441 PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
442 {
443   KSP_CG         *cg = (KSP_CG*)ksp->data;
444   PetscErrorCode ierr;
445   PetscBool      iascii;
446 
447   PetscFunctionBegin;
448   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
449   if (iascii) {
450 #if defined(PETSC_USE_COMPLEX)
451     ierr = PetscViewerASCIIPrintf(viewer,"  variant %s\n",KSPCGTypes[cg->type]);CHKERRQ(ierr);
452 #endif
453     if (cg->singlereduction) {
454       ierr = PetscViewerASCIIPrintf(viewer,"  using single-reduction variant\n");CHKERRQ(ierr);
455     }
456   }
457   PetscFunctionReturn(0);
458 }
459 
460 /*
461     KSPSetFromOptions_CG - Checks the options database for options related to the
462                            conjugate gradient method.
463 */
464 PetscErrorCode KSPSetFromOptions_CG(PetscOptionItems *PetscOptionsObject,KSP ksp)
465 {
466   PetscErrorCode ierr;
467   KSP_CG         *cg = (KSP_CG*)ksp->data;
468 
469   PetscFunctionBegin;
470   ierr = PetscOptionsHead(PetscOptionsObject,"KSP CG and CGNE options");CHKERRQ(ierr);
471 #if defined(PETSC_USE_COMPLEX)
472   ierr = PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
473                           (PetscEnum*)&cg->type,NULL);CHKERRQ(ierr);
474 #endif
475   ierr = PetscOptionsBool("-ksp_cg_single_reduction","Merge inner products into single MPIU_Allreduce()","KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,NULL);CHKERRQ(ierr);
476   ierr = PetscOptionsTail();CHKERRQ(ierr);
477   PetscFunctionReturn(0);
478 }
479 
480 /*
481     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
482                       This routine is registered below in KSPCreate_CG() and called from the
483                       routine KSPCGSetType() (see the file cgtype.c).
484 */
485 PetscErrorCode  KSPCGSetType_CG(KSP ksp,KSPCGType type)
486 {
487   KSP_CG *cg = (KSP_CG*)ksp->data;
488 
489   PetscFunctionBegin;
490   cg->type = type;
491   PetscFunctionReturn(0);
492 }
493 
494 /*
495     KSPCGUseSingleReduction_CG
496 
497     This routine sets a flag to use a variant of CG. Note that (in somewhat
498     atypical fashion) it also swaps out the routine called when KSPSolve()
499     is invoked.
500 */
501 static PetscErrorCode  KSPCGUseSingleReduction_CG(KSP ksp,PetscBool flg)
502 {
503   KSP_CG *cg = (KSP_CG*)ksp->data;
504 
505   PetscFunctionBegin;
506   cg->singlereduction = flg;
507   if (cg->singlereduction) {
508     ksp->ops->solve = KSPSolve_CG_SingleReduction;
509   } else {
510     ksp->ops->solve = KSPSolve_CG;
511   }
512   PetscFunctionReturn(0);
513 }
514 
515 /*
516     KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
517        function pointers for all the routines it needs to call (KSPSolve_CG() etc)
518 
519     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
520 */
521 /*MC
522      KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method
523 
524    Options Database Keys:
525 +   -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian, see KSPCGSetType()
526 .   -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric
527 -   -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single MPIU_Allreduce() call, see KSPCGUseSingleReduction()
528 
529    Level: beginner
530 
531    Notes:
532     The PCG method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite.
533 
534    Only left preconditioning is supported; there are several ways to motivate preconditioned CG, but they all produce the same algorithm.
535    One can interpret preconditioning A with B to mean any of the following\:
536 .n  (1) Solve a left-preconditioned system BAx = Bb, using inv(B) to define an inner product in the algorithm.
537 .n  (2) Solve a right-preconditioned system ABy = b, x = By, using B to define an inner product in the algorithm.
538 .n  (3) Solve a symmetrically-preconditioned system, E^TAEy = E^Tb, x = Ey, where B = EE^T.
539 .n  (4) Solve Ax=b with CG, but use the inner product defined by B to define the method [2].
540 .n  In all cases, the resulting algorithm only requires application of B to vectors.
541 
542    For complex numbers there are two different CG methods, one for Hermitian symmetric matrices and one for non-Hermitian symmetric matrices. Use
543    KSPCGSetType() to indicate which type you are using.
544 
545    Developer Notes:
546     KSPSolve_CG() should actually query the matrix to determine if it is Hermitian symmetric or not and NOT require the user to
547    indicate it to the KSP object.
548 
549    References:
550 .   1. - Magnus R. Hestenes and Eduard Stiefel, Methods of Conjugate Gradients for Solving Linear Systems,
551    Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
552 .   2. - Josef Malek and Zdenek Strakos, Preconditioning and the Conjugate Gradient Method in the Context of Solving PDEs,
553     SIAM, 2014.
554 
555 .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
556            KSPCGSetType(), KSPCGUseSingleReduction(), KSPPIPECG, KSPGROPPCG
557 
558 M*/
559 PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
560 {
561   PetscErrorCode ierr;
562   KSP_CG         *cg;
563 
564   PetscFunctionBegin;
565   ierr = PetscNewLog(ksp,&cg);CHKERRQ(ierr);
566 #if !defined(PETSC_USE_COMPLEX)
567   cg->type = KSP_CG_SYMMETRIC;
568 #else
569   cg->type = KSP_CG_HERMITIAN;
570 #endif
571   ksp->data = (void*)cg;
572 
573   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr);
574   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);CHKERRQ(ierr);
575   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);CHKERRQ(ierr);
576   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);CHKERRQ(ierr);
577 
578   /*
579        Sets the functions that are associated with this data structure
580        (in C++ this is the same as defining virtual functions)
581   */
582   ksp->ops->setup          = KSPSetUp_CG;
583   ksp->ops->solve          = KSPSolve_CG;
584   ksp->ops->destroy        = KSPDestroy_CG;
585   ksp->ops->view           = KSPView_CG;
586   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
587   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
588   ksp->ops->buildresidual  = KSPBuildResidualDefault;
589 
590   /*
591       Attach the function KSPCGSetType_CG() to this object. The routine
592       KSPCGSetType() checks for this attached function and calls it if it finds
593       it. (Sort of like a dynamic member function that can be added at run time
594   */
595   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",KSPCGSetType_CG);CHKERRQ(ierr);
596   ierr = PetscObjectComposeFunction((PetscObject)ksp,"KSPCGUseSingleReduction_C",KSPCGUseSingleReduction_CG);CHKERRQ(ierr);
597   PetscFunctionReturn(0);
598 }
599