xref: /petsc/src/ksp/ksp/impls/cg/cg.c (revision ae1ee55146a7ad071171b861759b85940c7e5c67)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith     This file implements the conjugate gradient method in PETSc as part of
34b9ad928SBarry Smith     KSP. You can use this as a starting point for implementing your own
44b9ad928SBarry Smith     Krylov method that is not provided with PETSc.
54b9ad928SBarry Smith 
64b9ad928SBarry Smith     The following basic routines are required for each Krylov method.
74b9ad928SBarry Smith         KSPCreate_XXX()          - Creates the Krylov context
84b9ad928SBarry Smith         KSPSetFromOptions_XXX()  - Sets runtime options
94b9ad928SBarry Smith         KSPSolve_XXX()           - Runs the Krylov method
104b9ad928SBarry Smith         KSPDestroy_XXX()         - Destroys the Krylov context, freeing all
114b9ad928SBarry Smith                                    memory it needed
124b9ad928SBarry Smith     Here the "_XXX" denotes a particular implementation, in this case
133cbce935SPierre Jolivet     we use _CG (e.g. KSPCreate_CG, KSPDestroy_CG). These routines
1463017de6SPatrick Sanan     are actually called via the common user interface routines
154b9ad928SBarry Smith     KSPSetType(), KSPSetFromOptions(), KSPSolve(), and KSPDestroy() so the
164b9ad928SBarry Smith     application code interface remains identical for all preconditioners.
174b9ad928SBarry Smith 
184b9ad928SBarry Smith     Other basic routines for the KSP objects include
194b9ad928SBarry Smith         KSPSetUp_XXX()
204b9ad928SBarry Smith         KSPView_XXX()            - Prints details of solver being used.
214b9ad928SBarry Smith 
22e94cfbe0SPatrick Sanan     Detailed Notes:
234b9ad928SBarry Smith     By default, this code implements the CG (Conjugate Gradient) method,
244b9ad928SBarry Smith     which is valid for real symmetric (and complex Hermitian) positive
254b9ad928SBarry Smith     definite matrices. Note that for the complex Hermitian case, the
264b9ad928SBarry Smith     VecDot() arguments within the code MUST remain in the order given
274b9ad928SBarry Smith     for correct computation of inner products.
284b9ad928SBarry Smith 
294b9ad928SBarry Smith     Reference: Hestenes and Steifel, 1952.
304b9ad928SBarry Smith 
314b9ad928SBarry Smith     By switching to the indefinite vector inner product, VecTDot(), the
324b9ad928SBarry Smith     same code is used for the complex symmetric case as well.  The user
334b9ad928SBarry Smith     must call KSPCGSetType(ksp,KSP_CG_SYMMETRIC) or use the option
341a952b2aSBarry Smith     -ksp_cg_type symmetric to invoke this variant for the complex case.
354b9ad928SBarry Smith     Note, however, that the complex symmetric code is NOT valid for
364b9ad928SBarry Smith     all such matrices ... and thus we don't recommend using this method.
374b9ad928SBarry Smith */
384b9ad928SBarry Smith /*
391f70db3aSBarry Smith     cgimpl.h defines the simple data structured used to store information
404b9ad928SBarry Smith     related to the type of matrix (e.g. complex symmetric) being solved and
410f9b5dd8SPierre Jolivet     data used during the optional Lanczos process used to compute eigenvalues
424b9ad928SBarry Smith */
43c6db04a5SJed Brown #include <../src/ksp/ksp/impls/cg/cgimpl.h> /*I "petscksp.h" I*/
4409573ac7SBarry Smith extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP, PetscReal *, PetscReal *);
4509573ac7SBarry Smith extern PetscErrorCode KSPComputeEigenvalues_CG(KSP, PetscInt, PetscReal *, PetscReal *, PetscInt *);
464b9ad928SBarry Smith 
KSPCGSetObjectiveTarget_CG(KSP ksp,PetscReal obj_min)47fb01098fSStefano Zampini static PetscErrorCode KSPCGSetObjectiveTarget_CG(KSP ksp, PetscReal obj_min)
48fb01098fSStefano Zampini {
49fb01098fSStefano Zampini   KSP_CG *cg = (KSP_CG *)ksp->data;
50fb01098fSStefano Zampini 
51fb01098fSStefano Zampini   PetscFunctionBegin;
52fb01098fSStefano Zampini   cg->obj_min = obj_min;
53fb01098fSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
54fb01098fSStefano Zampini }
55fb01098fSStefano Zampini 
KSPCGSetRadius_CG(KSP ksp,PetscReal radius)561a6f6928SStefano Zampini static PetscErrorCode KSPCGSetRadius_CG(KSP ksp, PetscReal radius)
571a6f6928SStefano Zampini {
581a6f6928SStefano Zampini   KSP_CG *cg = (KSP_CG *)ksp->data;
591a6f6928SStefano Zampini 
601a6f6928SStefano Zampini   PetscFunctionBegin;
611a6f6928SStefano Zampini   cg->radius = radius;
621a6f6928SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
631a6f6928SStefano Zampini }
641a6f6928SStefano Zampini 
KSPCGGetObjFcn_CG(KSP ksp,PetscReal * obj)65df8b9251SStefano Zampini static PetscErrorCode KSPCGGetObjFcn_CG(KSP ksp, PetscReal *obj)
66df8b9251SStefano Zampini {
67df8b9251SStefano Zampini   KSP_CG *cg = (KSP_CG *)ksp->data;
68df8b9251SStefano Zampini 
69df8b9251SStefano Zampini   PetscFunctionBegin;
70df8b9251SStefano Zampini   *obj = cg->obj;
71df8b9251SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
72df8b9251SStefano Zampini }
73df8b9251SStefano Zampini 
744b9ad928SBarry Smith /*
754b9ad928SBarry Smith      KSPSetUp_CG - Sets up the workspace needed by the CG method.
764b9ad928SBarry Smith 
774b9ad928SBarry Smith       This is called once, usually automatically by KSPSolve() or KSPSetUp()
784b9ad928SBarry Smith      but can be called directly by KSPSetUp()
794b9ad928SBarry Smith */
KSPSetUp_CG(KSP ksp)80d71ae5a4SJacob Faibussowitsch static PetscErrorCode KSPSetUp_CG(KSP ksp)
81d71ae5a4SJacob Faibussowitsch {
824b9ad928SBarry Smith   KSP_CG  *cgP   = (KSP_CG *)ksp->data;
83608df7e2SBarry Smith   PetscInt maxit = ksp->max_it, nwork = 3;
844b9ad928SBarry Smith 
854b9ad928SBarry Smith   PetscFunctionBegin;
864b9ad928SBarry Smith   /* get work vectors needed by CG */
87043bf087SBarry Smith   if (cgP->singlereduction) nwork += 2;
889566063dSJacob Faibussowitsch   PetscCall(KSPSetWorkVecs(ksp, nwork));
894b9ad928SBarry Smith 
904b9ad928SBarry Smith   /*
914802679fSStefano Zampini      If user requested computations of eigenvalues then allocate
924b9ad928SBarry Smith      work space needed
934b9ad928SBarry Smith   */
944b9ad928SBarry Smith   if (ksp->calc_sings) {
959566063dSJacob Faibussowitsch     PetscCall(PetscFree4(cgP->e, cgP->d, cgP->ee, cgP->dd));
96a25662f8SPablo Brubeck     PetscCall(PetscMalloc4(maxit + 1, &cgP->e, maxit, &cgP->d, maxit, &cgP->ee, maxit, &cgP->dd));
97b30eab52SKarl Rupp 
984b9ad928SBarry Smith     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
994b9ad928SBarry Smith     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
1004b9ad928SBarry Smith   }
1013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1024b9ad928SBarry Smith }
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith /*
1058d36b627SPatrick Sanan      A macro used in the following KSPSolve_CG and KSPSolve_CG_SingleReduction routines
10614579225SPatrick Sanan */
1074ad8454bSPierre Jolivet #define VecXDot(x, y, a) (cg->type == KSP_CG_HERMITIAN ? VecDot(x, y, a) : VecTDot(x, y, a))
10814579225SPatrick Sanan 
10914579225SPatrick Sanan /*
11006e420deSBarry Smith      KSPSolve_CG - This routine actually applies the conjugate gradient method
11106e420deSBarry Smith 
11214579225SPatrick Sanan      Note : this routine can be replaced with another one (see below) which implements
11314579225SPatrick Sanan             another variant of CG.
1144b9ad928SBarry Smith 
1154b9ad928SBarry Smith    Input Parameter:
1164b9ad928SBarry Smith .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
1174b9ad928SBarry Smith             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
1184b9ad928SBarry Smith */
KSPSolve_CG(KSP ksp)119d71ae5a4SJacob Faibussowitsch static PetscErrorCode KSPSolve_CG(KSP ksp)
120d71ae5a4SJacob Faibussowitsch {
121e10ebe35SBarry Smith   PetscInt    i, stored_max_it, eigs;
1220a545947SLisandro Dalcin   PetscScalar dpi = 0.0, a = 1.0, beta, betaold = 1.0, b = 0, *e = NULL, *d = NULL, dpiold;
1234b9ad928SBarry Smith   PetscReal   dp = 0.0;
1241a6f6928SStefano Zampini   PetscReal   r2, norm_p, norm_d, dMp;
12563017de6SPatrick Sanan   Vec         X, B, Z, R, P, W;
1264b9ad928SBarry Smith   KSP_CG     *cg;
1274b9ad928SBarry Smith   Mat         Amat, Pmat;
128fb01098fSStefano Zampini   PetscBool   diagonalscale, testobj;
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   PetscFunctionBegin;
1319566063dSJacob Faibussowitsch   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
1325f80ce2aSJacob Faibussowitsch   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   cg            = (KSP_CG *)ksp->data;
1354b9ad928SBarry Smith   eigs          = ksp->calc_sings;
1364b9ad928SBarry Smith   stored_max_it = ksp->max_it;
1374b9ad928SBarry Smith   X             = ksp->vec_sol;
1384b9ad928SBarry Smith   B             = ksp->vec_rhs;
1394b9ad928SBarry Smith   R             = ksp->work[0];
1404b9ad928SBarry Smith   Z             = ksp->work[1];
1414b9ad928SBarry Smith   P             = ksp->work[2];
142043bf087SBarry Smith   W             = Z;
1431a6f6928SStefano Zampini   r2            = PetscSqr(cg->radius);
1444b9ad928SBarry Smith 
1459371c9d4SSatish Balay   if (eigs) {
1469371c9d4SSatish Balay     e    = cg->e;
1479371c9d4SSatish Balay     d    = cg->d;
1489371c9d4SSatish Balay     e[0] = 0.0;
1499371c9d4SSatish Balay   }
1509566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith   ksp->its = 0;
1534b9ad928SBarry Smith   if (!ksp->guess_zero) {
1549566063dSJacob Faibussowitsch     PetscCall(KSP_MatMult(ksp, Amat, X, R)); /*    r <- b - Ax                       */
155fb01098fSStefano Zampini 
1569566063dSJacob Faibussowitsch     PetscCall(VecAYPX(R, -1.0, B));
157fb01098fSStefano Zampini     if (cg->radius) { /* XXX direction? */
1581a6f6928SStefano Zampini       PetscCall(VecNorm(X, NORM_2, &norm_d));
1591a6f6928SStefano Zampini       norm_d *= norm_d;
1601a6f6928SStefano Zampini     }
1614b9ad928SBarry Smith   } else {
1629566063dSJacob Faibussowitsch     PetscCall(VecCopy(B, R)); /*    r <- b (x is 0)                   */
1631a6f6928SStefano Zampini     norm_d = 0.0;
1644b9ad928SBarry Smith   }
165d8256aecSBarry Smith   /* This may be true only on a subset of MPI ranks; setting it here so it will be detected by the first norm computation below */
166f480ea8aSBarry Smith   PetscCall(VecFlag(R, ksp->reason == KSP_DIVERGED_PC_FAILED));
167d33d232eStmunson 
1689e568555SJed Brown   switch (ksp->normtype) {
1699e568555SJed Brown   case KSP_NORM_PRECONDITIONED:
1709566063dSJacob Faibussowitsch     PetscCall(KSP_PCApply(ksp, R, Z));  /*    z <- Br                           */
1719566063dSJacob Faibussowitsch     PetscCall(VecNorm(Z, NORM_2, &dp)); /*    dp <- z'*z = e'*A'*B'*B*A*e       */
172c0decd05SBarry Smith     KSPCheckNorm(ksp, dp);
1739e568555SJed Brown     break;
1749e568555SJed Brown   case KSP_NORM_UNPRECONDITIONED:
1759566063dSJacob Faibussowitsch     PetscCall(VecNorm(R, NORM_2, &dp)); /*    dp <- r'*r = e'*A'*A*e            */
176c0decd05SBarry Smith     KSPCheckNorm(ksp, dp);
1779e568555SJed Brown     break;
1789e568555SJed Brown   case KSP_NORM_NATURAL:
1799566063dSJacob Faibussowitsch     PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
1809566063dSJacob Faibussowitsch     PetscCall(VecXDot(Z, R, &beta));   /*    beta <- z'*r                      */
181422a814eSBarry Smith     KSPCheckDot(ksp, beta);
1828f1a2a5eSBarry Smith     dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
1839e568555SJed Brown     break;
184d71ae5a4SJacob Faibussowitsch   case KSP_NORM_NONE:
185d71ae5a4SJacob Faibussowitsch     dp = 0.0;
186d71ae5a4SJacob Faibussowitsch     break;
187d71ae5a4SJacob Faibussowitsch   default:
188d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
1899e568555SJed Brown   }
190fb01098fSStefano Zampini 
191fb01098fSStefano Zampini   /* Initialize objective function
192fb01098fSStefano Zampini      obj = 1/2 x^T A x - x^T b */
193fb01098fSStefano Zampini   testobj = (PetscBool)(cg->obj_min < 0.0);
194fb01098fSStefano Zampini   PetscCall(VecXDot(R, X, &a));
195fb01098fSStefano Zampini   cg->obj = 0.5 * PetscRealPart(a);
196fb01098fSStefano Zampini   PetscCall(VecXDot(B, X, &a));
197fb01098fSStefano Zampini   cg->obj -= 0.5 * PetscRealPart(a);
198fb01098fSStefano Zampini 
199fc69cb9eSStefano Zampini   if (testobj) PetscCall(PetscInfo(ksp, "it %" PetscInt_FMT " obj %g\n", ksp->its, (double)cg->obj));
2009566063dSJacob Faibussowitsch   PetscCall(KSPLogResidualHistory(ksp, dp));
201fb01098fSStefano Zampini   PetscCall(KSPMonitor(ksp, ksp->its, dp));
2024b9ad928SBarry Smith   ksp->rnorm = dp;
2034b9ad928SBarry Smith 
204fb01098fSStefano Zampini   PetscCall((*ksp->converged)(ksp, ksp->its, dp, &ksp->reason, ksp->cnvP)); /* test for convergence */
205fb01098fSStefano Zampini 
206fb01098fSStefano Zampini   if (!ksp->reason && testobj && cg->obj <= cg->obj_min) {
207fb01098fSStefano Zampini     PetscCall(PetscInfo(ksp, "converged to objective target minimum\n"));
208fb01098fSStefano Zampini     ksp->reason = KSP_CONVERGED_ATOL;
209fb01098fSStefano Zampini   }
210fb01098fSStefano Zampini 
2113ba16761SJacob Faibussowitsch   if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS);
2123b020e0cSBarry Smith 
213789736e1SBarry Smith   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) PetscCall(KSP_PCApply(ksp, R, Z)); /*     z <- Br                           */
214fc165ae2SBarry Smith   if (ksp->normtype != KSP_NORM_NATURAL) {
2159566063dSJacob Faibussowitsch     PetscCall(VecXDot(Z, R, &beta)); /*     beta <- z'*r                      */
216422a814eSBarry Smith     KSPCheckDot(ksp, beta);
217fc165ae2SBarry Smith   }
218fc165ae2SBarry Smith 
2194b9ad928SBarry Smith   i = 0;
2204b9ad928SBarry Smith   do {
2214b9ad928SBarry Smith     ksp->its = i + 1;
2224b9ad928SBarry Smith     if (beta == 0.0) {
2234b9ad928SBarry Smith       ksp->reason = KSP_CONVERGED_ATOL;
2249566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ksp, "converged due to beta = 0\n"));
2254b9ad928SBarry Smith       break;
2264b9ad928SBarry Smith #if !defined(PETSC_USE_COMPLEX)
227e0a67c7fSBarry Smith     } else if ((i > 0) && (beta * betaold < 0.0)) {
228a25662f8SPablo Brubeck       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite preconditioner, beta %g, betaold %g", (double)PetscRealPart(beta), (double)PetscRealPart(betaold));
2294b9ad928SBarry Smith       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
2309566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ksp, "diverging due to indefinite preconditioner\n"));
2314b9ad928SBarry Smith       break;
2324b9ad928SBarry Smith #endif
2334b9ad928SBarry Smith     }
2344b9ad928SBarry Smith     if (!i) {
2359566063dSJacob Faibussowitsch       PetscCall(VecCopy(Z, P)); /*     p <- z                           */
2361a6f6928SStefano Zampini       if (cg->radius) {
2371a6f6928SStefano Zampini         PetscCall(VecNorm(P, NORM_2, &norm_p));
2381a6f6928SStefano Zampini         norm_p *= norm_p;
2391a6f6928SStefano Zampini         dMp = 0.0;
2403a7d0413SPierre Jolivet         if (!ksp->guess_zero) PetscCall(VecDotRealPart(X, P, &dMp));
2411a6f6928SStefano Zampini       }
2424869883bSvictor       b = 0.0;
2434b9ad928SBarry Smith     } else {
2444b9ad928SBarry Smith       b = beta / betaold;
2454b9ad928SBarry Smith       if (eigs) {
2465f80ce2aSJacob Faibussowitsch         PetscCheck(ksp->max_it == stored_max_it, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Can not change maxit AND calculate eigenvalues");
2478f1a2a5eSBarry Smith         e[i] = PetscSqrtReal(PetscAbsScalar(b)) / a;
2484b9ad928SBarry Smith       }
2499566063dSJacob Faibussowitsch       PetscCall(VecAYPX(P, b, Z)); /*     p <- z + b* p                    */
2501a6f6928SStefano Zampini       if (cg->radius) {
2511a6f6928SStefano Zampini         PetscCall(VecDotRealPart(X, P, &dMp));
2521a6f6928SStefano Zampini         PetscCall(VecNorm(P, NORM_2, &norm_p));
2531a6f6928SStefano Zampini         norm_p *= norm_p;
2541a6f6928SStefano Zampini       }
2554b9ad928SBarry Smith     }
256608df7e2SBarry Smith     dpiold = dpi;
2579566063dSJacob Faibussowitsch     PetscCall(KSP_MatMult(ksp, Amat, P, W)); /*     w <- Ap                          */
2589566063dSJacob Faibussowitsch     PetscCall(VecXDot(P, W, &dpi));          /*     dpi <- p'w                       */
259647c6a20SBarry Smith     KSPCheckDot(ksp, dpi);
26063017de6SPatrick Sanan     betaold = beta;
26163017de6SPatrick Sanan 
262c0decd05SBarry Smith     if ((dpi == 0.0) || ((i > 0) && ((PetscSign(PetscRealPart(dpi)) * PetscSign(PetscRealPart(dpiold))) < 0.0))) {
2631a6f6928SStefano Zampini       if (cg->radius) {
2641a6f6928SStefano Zampini         a = 0.0;
2651a6f6928SStefano Zampini         if (i == 0) {
2661a6f6928SStefano Zampini           if (norm_p > 0.0) {
2671a6f6928SStefano Zampini             a = PetscSqrtReal(r2 / norm_p);
2681a6f6928SStefano Zampini           } else {
2691a6f6928SStefano Zampini             PetscCall(VecNorm(R, NORM_2, &dp));
2701a6f6928SStefano Zampini             a = cg->radius > dp ? 1.0 : cg->radius / dp;
2711a6f6928SStefano Zampini           }
2721a6f6928SStefano Zampini         } else if (norm_p > 0.0) {
2731a6f6928SStefano Zampini           a = (PetscSqrtReal(dMp * dMp + norm_p * (r2 - norm_d)) - dMp) / norm_p;
2741a6f6928SStefano Zampini         }
2751a6f6928SStefano Zampini         PetscCall(VecAXPY(X, a, P)); /*     x <- x + ap                      */
276fb01098fSStefano Zampini         cg->obj += PetscRealPart(a * (0.5 * a * dpi - betaold));
2772a28d964SStefano Zampini       }
278fc69cb9eSStefano Zampini       if (testobj) PetscCall(PetscInfo(ksp, "it %" PetscInt_FMT " N obj %g\n", i + 1, (double)cg->obj));
2792a28d964SStefano Zampini       if (ksp->converged_neg_curve) {
2801a6f6928SStefano Zampini         PetscCall(PetscInfo(ksp, "converged due to negative curvature: %g\n", (double)(PetscRealPart(dpi))));
2811a6f6928SStefano Zampini         ksp->reason = KSP_CONVERGED_NEG_CURVE;
2821a6f6928SStefano Zampini       } else {
2835f80ce2aSJacob Faibussowitsch         PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite matrix, dpi %g, dpiold %g", (double)PetscRealPart(dpi), (double)PetscRealPart(dpiold));
28463017de6SPatrick Sanan         ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
2851a6f6928SStefano Zampini         PetscCall(PetscInfo(ksp, "diverging due to indefinite matrix\n"));
2861a6f6928SStefano Zampini       }
28763017de6SPatrick Sanan       break;
28863017de6SPatrick Sanan     }
28963017de6SPatrick Sanan     a = beta / dpi; /*     a = beta/p'w                     */
29063017de6SPatrick Sanan     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b)) * e[i] + 1.0 / a;
291fb01098fSStefano Zampini     if (cg->radius) { /* Steihaugh-Toint */
2921a6f6928SStefano Zampini       PetscReal norm_dp1 = norm_d + PetscRealPart(a) * (2.0 * dMp + PetscRealPart(a) * norm_p);
2931a6f6928SStefano Zampini       if (norm_dp1 > r2) {
2941a6f6928SStefano Zampini         ksp->reason = KSP_CONVERGED_STEP_LENGTH;
2951a6f6928SStefano Zampini         PetscCall(PetscInfo(ksp, "converged to the trust region radius %g\n", (double)cg->radius));
2961a6f6928SStefano Zampini         if (norm_p > 0.0) {
2971a6f6928SStefano Zampini           dp = (PetscSqrtReal(dMp * dMp + norm_p * (r2 - norm_d)) - dMp) / norm_p;
2981a6f6928SStefano Zampini           PetscCall(VecAXPY(X, dp, P)); /*     x <- x + ap                      */
299fb01098fSStefano Zampini           cg->obj += PetscRealPart(dp * (0.5 * dp * dpi - beta));
3001a6f6928SStefano Zampini         }
301fc69cb9eSStefano Zampini         if (testobj) PetscCall(PetscInfo(ksp, "it %" PetscInt_FMT " R obj %g\n", i + 1, (double)cg->obj));
3021a6f6928SStefano Zampini         break;
3031a6f6928SStefano Zampini       }
3041a6f6928SStefano Zampini     }
3059566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X, a, P));  /*     x <- x + ap                      */
3069566063dSJacob Faibussowitsch     PetscCall(VecAXPY(R, -a, W)); /*     r <- r - aw                      */
30763017de6SPatrick Sanan     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i + 2) {
3089566063dSJacob Faibussowitsch       PetscCall(KSP_PCApply(ksp, R, Z));  /*     z <- Br                          */
3099566063dSJacob Faibussowitsch       PetscCall(VecNorm(Z, NORM_2, &dp)); /*     dp <- z'*z                       */
310c0decd05SBarry Smith       KSPCheckNorm(ksp, dp);
31163017de6SPatrick Sanan     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i + 2) {
3129566063dSJacob Faibussowitsch       PetscCall(VecNorm(R, NORM_2, &dp)); /*     dp <- r'*r                       */
313c0decd05SBarry Smith       KSPCheckNorm(ksp, dp);
31463017de6SPatrick Sanan     } else if (ksp->normtype == KSP_NORM_NATURAL) {
3159566063dSJacob Faibussowitsch       PetscCall(KSP_PCApply(ksp, R, Z)); /*     z <- Br                          */
3169566063dSJacob Faibussowitsch       PetscCall(VecXDot(Z, R, &beta));   /*     beta <- r'*z                     */
31763017de6SPatrick Sanan       KSPCheckDot(ksp, beta);
31863017de6SPatrick Sanan       dp = PetscSqrtReal(PetscAbsScalar(beta));
31963017de6SPatrick Sanan     } else {
32063017de6SPatrick Sanan       dp = 0.0;
32163017de6SPatrick Sanan     }
322fb01098fSStefano Zampini     cg->obj -= PetscRealPart(0.5 * a * betaold);
323fc69cb9eSStefano Zampini     if (testobj) PetscCall(PetscInfo(ksp, "it %" PetscInt_FMT " obj %g\n", i + 1, (double)cg->obj));
324fb01098fSStefano Zampini 
32563017de6SPatrick Sanan     ksp->rnorm = dp;
3269566063dSJacob Faibussowitsch     PetscCall(KSPLogResidualHistory(ksp, dp));
3279566063dSJacob Faibussowitsch     PetscCall(KSPMonitor(ksp, i + 1, dp));
3289566063dSJacob Faibussowitsch     PetscCall((*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP));
329fb01098fSStefano Zampini 
330fb01098fSStefano Zampini     if (!ksp->reason && testobj && cg->obj <= cg->obj_min) {
331fb01098fSStefano Zampini       PetscCall(PetscInfo(ksp, "converged to objective target minimum\n"));
332fb01098fSStefano Zampini       ksp->reason = KSP_CONVERGED_ATOL;
333fb01098fSStefano Zampini     }
334fb01098fSStefano Zampini 
33563017de6SPatrick Sanan     if (ksp->reason) break;
33663017de6SPatrick Sanan 
3371a6f6928SStefano Zampini     if (cg->radius) {
3381a6f6928SStefano Zampini       PetscCall(VecNorm(X, NORM_2, &norm_d));
3391a6f6928SStefano Zampini       norm_d *= norm_d;
3401a6f6928SStefano Zampini     }
3411a6f6928SStefano Zampini 
342b6555650SPierre Jolivet     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && ksp->normtype != KSP_NORM_NATURAL) || ksp->chknorm >= i + 2) PetscCall(KSP_PCApply(ksp, R, Z)); /*     z <- Br                          */
343b6555650SPierre Jolivet     if (ksp->normtype != KSP_NORM_NATURAL || ksp->chknorm >= i + 2) {
3449566063dSJacob Faibussowitsch       PetscCall(VecXDot(Z, R, &beta)); /*     beta <- z'*r                     */
34563017de6SPatrick Sanan       KSPCheckDot(ksp, beta);
34663017de6SPatrick Sanan     }
34763017de6SPatrick Sanan 
34863017de6SPatrick Sanan     i++;
34963017de6SPatrick Sanan   } while (i < ksp->max_it);
35063017de6SPatrick Sanan   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35263017de6SPatrick Sanan }
35363017de6SPatrick Sanan 
35463017de6SPatrick Sanan /*
35563017de6SPatrick Sanan        KSPSolve_CG_SingleReduction
35663017de6SPatrick Sanan 
35763017de6SPatrick Sanan        This variant of CG is identical in exact arithmetic to the standard algorithm,
35863017de6SPatrick Sanan        but is rearranged to use only a single reduction stage per iteration, using additional
35963017de6SPatrick Sanan        intermediate vectors.
36014579225SPatrick Sanan 
36114579225SPatrick Sanan        See KSPCGUseSingleReduction_CG()
36214579225SPatrick Sanan 
36363017de6SPatrick Sanan */
KSPSolve_CG_SingleReduction(KSP ksp)364d71ae5a4SJacob Faibussowitsch static PetscErrorCode KSPSolve_CG_SingleReduction(KSP ksp)
365d71ae5a4SJacob Faibussowitsch {
36663017de6SPatrick Sanan   PetscInt    i, stored_max_it, eigs;
3670a545947SLisandro Dalcin   PetscScalar dpi = 0.0, a = 1.0, beta, betaold = 1.0, b = 0, *e = NULL, *d = NULL, delta, dpiold, tmp[2];
36863017de6SPatrick Sanan   PetscReal   dp = 0.0;
36963017de6SPatrick Sanan   Vec         X, B, Z, R, P, S, W, tmpvecs[2];
37063017de6SPatrick Sanan   KSP_CG     *cg;
37163017de6SPatrick Sanan   Mat         Amat, Pmat;
37263017de6SPatrick Sanan   PetscBool   diagonalscale;
37363017de6SPatrick Sanan 
37463017de6SPatrick Sanan   PetscFunctionBegin;
3759566063dSJacob Faibussowitsch   PetscCall(PCGetDiagonalScale(ksp->pc, &diagonalscale));
3765f80ce2aSJacob Faibussowitsch   PetscCheck(!diagonalscale, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Krylov method %s does not support diagonal scaling", ((PetscObject)ksp)->type_name);
37763017de6SPatrick Sanan 
37863017de6SPatrick Sanan   cg            = (KSP_CG *)ksp->data;
37963017de6SPatrick Sanan   eigs          = ksp->calc_sings;
38063017de6SPatrick Sanan   stored_max_it = ksp->max_it;
38163017de6SPatrick Sanan   X             = ksp->vec_sol;
38263017de6SPatrick Sanan   B             = ksp->vec_rhs;
38363017de6SPatrick Sanan   R             = ksp->work[0];
38463017de6SPatrick Sanan   Z             = ksp->work[1];
38563017de6SPatrick Sanan   P             = ksp->work[2];
38663017de6SPatrick Sanan   S             = ksp->work[3];
38763017de6SPatrick Sanan   W             = ksp->work[4];
38863017de6SPatrick Sanan 
3899371c9d4SSatish Balay   if (eigs) {
3909371c9d4SSatish Balay     e    = cg->e;
3919371c9d4SSatish Balay     d    = cg->d;
3929371c9d4SSatish Balay     e[0] = 0.0;
3939371c9d4SSatish Balay   }
3949566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(ksp->pc, &Amat, &Pmat));
39563017de6SPatrick Sanan 
39663017de6SPatrick Sanan   ksp->its = 0;
39763017de6SPatrick Sanan   if (!ksp->guess_zero) {
3989566063dSJacob Faibussowitsch     PetscCall(KSP_MatMult(ksp, Amat, X, R)); /*    r <- b - Ax                       */
3999566063dSJacob Faibussowitsch     PetscCall(VecAYPX(R, -1.0, B));
40063017de6SPatrick Sanan   } else {
4019566063dSJacob Faibussowitsch     PetscCall(VecCopy(B, R)); /*    r <- b (x is 0)                   */
40263017de6SPatrick Sanan   }
40363017de6SPatrick Sanan 
40463017de6SPatrick Sanan   switch (ksp->normtype) {
40563017de6SPatrick Sanan   case KSP_NORM_PRECONDITIONED:
4069566063dSJacob Faibussowitsch     PetscCall(KSP_PCApply(ksp, R, Z));  /*    z <- Br                           */
4079566063dSJacob Faibussowitsch     PetscCall(VecNorm(Z, NORM_2, &dp)); /*    dp <- z'*z = e'*A'*B'*B*A'*e'     */
408c0decd05SBarry Smith     KSPCheckNorm(ksp, dp);
40963017de6SPatrick Sanan     break;
41063017de6SPatrick Sanan   case KSP_NORM_UNPRECONDITIONED:
4119566063dSJacob Faibussowitsch     PetscCall(VecNorm(R, NORM_2, &dp)); /*    dp <- r'*r = e'*A'*A*e            */
412c0decd05SBarry Smith     KSPCheckNorm(ksp, dp);
41363017de6SPatrick Sanan     break;
41463017de6SPatrick Sanan   case KSP_NORM_NATURAL:
4159566063dSJacob Faibussowitsch     PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
4169566063dSJacob Faibussowitsch     PetscCall(KSP_MatMult(ksp, Amat, Z, S));
4179566063dSJacob Faibussowitsch     PetscCall(VecXDot(Z, S, &delta)); /*    delta <- z'*A*z = r'*B*A*B*r      */
4189566063dSJacob Faibussowitsch     PetscCall(VecXDot(Z, R, &beta));  /*    beta <- z'*r                      */
41963017de6SPatrick Sanan     KSPCheckDot(ksp, beta);
42063017de6SPatrick Sanan     dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
42163017de6SPatrick Sanan     break;
422d71ae5a4SJacob Faibussowitsch   case KSP_NORM_NONE:
423d71ae5a4SJacob Faibussowitsch     dp = 0.0;
424d71ae5a4SJacob Faibussowitsch     break;
425d71ae5a4SJacob Faibussowitsch   default:
426d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "%s", KSPNormTypes[ksp->normtype]);
42763017de6SPatrick Sanan   }
4289566063dSJacob Faibussowitsch   PetscCall(KSPLogResidualHistory(ksp, dp));
4299566063dSJacob Faibussowitsch   PetscCall(KSPMonitor(ksp, 0, dp));
43063017de6SPatrick Sanan   ksp->rnorm = dp;
43163017de6SPatrick Sanan 
4329566063dSJacob Faibussowitsch   PetscCall((*ksp->converged)(ksp, 0, dp, &ksp->reason, ksp->cnvP)); /* test for convergence */
4333ba16761SJacob Faibussowitsch   if (ksp->reason) PetscFunctionReturn(PETSC_SUCCESS);
43463017de6SPatrick Sanan 
435ac530a7eSPierre Jolivet   if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
43663017de6SPatrick Sanan   if (ksp->normtype != KSP_NORM_NATURAL) {
4379566063dSJacob Faibussowitsch     PetscCall(KSP_MatMult(ksp, Amat, Z, S));
4389566063dSJacob Faibussowitsch     PetscCall(VecXDot(Z, S, &delta)); /*    delta <- z'*A*z = r'*B*A*B*r      */
4399566063dSJacob Faibussowitsch     PetscCall(VecXDot(Z, R, &beta));  /*    beta <- z'*r                      */
44063017de6SPatrick Sanan     KSPCheckDot(ksp, beta);
44163017de6SPatrick Sanan   }
44263017de6SPatrick Sanan 
44363017de6SPatrick Sanan   i = 0;
44463017de6SPatrick Sanan   do {
44563017de6SPatrick Sanan     ksp->its = i + 1;
44663017de6SPatrick Sanan     if (beta == 0.0) {
44763017de6SPatrick Sanan       ksp->reason = KSP_CONVERGED_ATOL;
4489566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ksp, "converged due to beta = 0\n"));
44963017de6SPatrick Sanan       break;
45063017de6SPatrick Sanan #if !defined(PETSC_USE_COMPLEX)
45163017de6SPatrick Sanan     } else if ((i > 0) && (beta * betaold < 0.0)) {
4525f80ce2aSJacob Faibussowitsch       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite preconditioner");
45363017de6SPatrick Sanan       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
4549566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ksp, "diverging due to indefinite preconditioner\n"));
45563017de6SPatrick Sanan       break;
45663017de6SPatrick Sanan #endif
45763017de6SPatrick Sanan     }
45863017de6SPatrick Sanan     if (!i) {
4599566063dSJacob Faibussowitsch       PetscCall(VecCopy(Z, P)); /*    p <- z                           */
46063017de6SPatrick Sanan       b = 0.0;
46163017de6SPatrick Sanan     } else {
46263017de6SPatrick Sanan       b = beta / betaold;
46363017de6SPatrick Sanan       if (eigs) {
4645f80ce2aSJacob Faibussowitsch         PetscCheck(ksp->max_it == stored_max_it, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "Can not change maxit AND calculate eigenvalues");
46563017de6SPatrick Sanan         e[i] = PetscSqrtReal(PetscAbsScalar(b)) / a;
46663017de6SPatrick Sanan       }
4679566063dSJacob Faibussowitsch       PetscCall(VecAYPX(P, b, Z)); /*    p <- z + b* p                     */
46863017de6SPatrick Sanan     }
46963017de6SPatrick Sanan     dpiold = dpi;
47063017de6SPatrick Sanan     if (!i) {
4719566063dSJacob Faibussowitsch       PetscCall(KSP_MatMult(ksp, Amat, P, W)); /*    w <- Ap                           */
4729566063dSJacob Faibussowitsch       PetscCall(VecXDot(P, W, &dpi));          /*    dpi <- p'w                        */
473043bf087SBarry Smith     } else {
4749566063dSJacob Faibussowitsch       PetscCall(VecAYPX(W, beta / betaold, S));                 /*    w <- Ap                           */
475d8b2d19eSBarry Smith       dpi = delta - beta * beta * dpiold / (betaold * betaold); /*    dpi <- p'w                        */
476043bf087SBarry Smith     }
477608df7e2SBarry Smith     betaold = beta;
478422a814eSBarry Smith     KSPCheckDot(ksp, beta);
479d33d232eStmunson 
480e0a67c7fSBarry Smith     if ((dpi == 0.0) || ((i > 0) && (PetscRealPart(dpi * dpiold) <= 0.0))) {
4815f80ce2aSJacob Faibussowitsch       PetscCheck(!ksp->errorifnotconverged, PetscObjectComm((PetscObject)ksp), PETSC_ERR_NOT_CONVERGED, "Diverged due to indefinite matrix");
4826aee1118SBarry Smith       ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
4839566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ksp, "diverging due to indefinite or negative definite matrix\n"));
4846aee1118SBarry Smith       break;
4856aee1118SBarry Smith     }
486043bf087SBarry Smith     a = beta / dpi; /*    a = beta/p'w                      */
487b30eab52SKarl Rupp     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b)) * e[i] + 1.0 / a;
4889566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X, a, P));  /*    x <- x + ap                       */
4899566063dSJacob Faibussowitsch     PetscCall(VecAXPY(R, -a, W)); /*    r <- r - aw                       */
490dc2b3c31SBarry Smith     if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i + 2) {
4919566063dSJacob Faibussowitsch       PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
4929566063dSJacob Faibussowitsch       PetscCall(KSP_MatMult(ksp, Amat, Z, S));
4939566063dSJacob Faibussowitsch       PetscCall(VecNorm(Z, NORM_2, &dp)); /*    dp <- z'*z                        */
494c0decd05SBarry Smith       KSPCheckNorm(ksp, dp);
495dc2b3c31SBarry Smith     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i + 2) {
4969566063dSJacob Faibussowitsch       PetscCall(VecNorm(R, NORM_2, &dp)); /*    dp <- r'*r                        */
497c0decd05SBarry Smith       KSPCheckNorm(ksp, dp);
498ce9499c7SBarry Smith     } else if (ksp->normtype == KSP_NORM_NATURAL) {
4999566063dSJacob Faibussowitsch       PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
5009371c9d4SSatish Balay       tmpvecs[0] = S;
5019371c9d4SSatish Balay       tmpvecs[1] = R;
5029566063dSJacob Faibussowitsch       PetscCall(KSP_MatMult(ksp, Amat, Z, S));
5039566063dSJacob Faibussowitsch       PetscCall(VecMDot(Z, 2, tmpvecs, tmp)); /*    delta <- z'*A*z = r'*B*A*B*r      */
5049371c9d4SSatish Balay       delta = tmp[0];
5059371c9d4SSatish Balay       beta  = tmp[1]; /*    beta <- z'*r                      */
506422a814eSBarry Smith       KSPCheckDot(ksp, beta);
50763017de6SPatrick Sanan       dp = PetscSqrtReal(PetscAbsScalar(beta)); /*    dp <- r'*z = r'*B*r = e'*A'*B*A*e */
5084b9ad928SBarry Smith     } else {
5094b9ad928SBarry Smith       dp = 0.0;
5104b9ad928SBarry Smith     }
5114b9ad928SBarry Smith     ksp->rnorm = dp;
5129566063dSJacob Faibussowitsch     PetscCall(KSPLogResidualHistory(ksp, dp));
5139566063dSJacob Faibussowitsch     PetscCall(KSPMonitor(ksp, i + 1, dp));
5149566063dSJacob Faibussowitsch     PetscCall((*ksp->converged)(ksp, i + 1, dp, &ksp->reason, ksp->cnvP));
5154b9ad928SBarry Smith     if (ksp->reason) break;
516fc165ae2SBarry Smith 
517b6555650SPierre Jolivet     if ((ksp->normtype != KSP_NORM_PRECONDITIONED && ksp->normtype != KSP_NORM_NATURAL) || ksp->chknorm >= i + 2) {
5189566063dSJacob Faibussowitsch       PetscCall(KSP_PCApply(ksp, R, Z)); /*    z <- Br                           */
5199566063dSJacob Faibussowitsch       PetscCall(KSP_MatMult(ksp, Amat, Z, S));
520608df7e2SBarry Smith     }
521b6555650SPierre Jolivet     if (ksp->normtype != KSP_NORM_NATURAL || ksp->chknorm >= i + 2) {
5229371c9d4SSatish Balay       tmpvecs[0] = S;
5239371c9d4SSatish Balay       tmpvecs[1] = R;
5249566063dSJacob Faibussowitsch       PetscCall(VecMDot(Z, 2, tmpvecs, tmp));
5259371c9d4SSatish Balay       delta = tmp[0];
5269371c9d4SSatish Balay       beta  = tmp[1];         /*    delta <- z'*A*z = r'*B'*A*B*r     */
52763017de6SPatrick Sanan       KSPCheckDot(ksp, beta); /*    beta <- z'*r                      */
528fc165ae2SBarry Smith     }
529fc165ae2SBarry Smith 
5304b9ad928SBarry Smith     i++;
5314b9ad928SBarry Smith   } while (i < ksp->max_it);
532b30eab52SKarl Rupp   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
5333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5344b9ad928SBarry Smith }
5354b9ad928SBarry Smith 
53663017de6SPatrick Sanan /*
53763017de6SPatrick Sanan      KSPDestroy_CG - Frees resources allocated in KSPSetup_CG and clears function
53863017de6SPatrick Sanan                      compositions from KSPCreate_CG. If adding your own KSP implementation,
53963017de6SPatrick Sanan                      you must be sure to free all allocated resources here to prevent
54063017de6SPatrick Sanan                      leaks.
54163017de6SPatrick Sanan */
KSPDestroy_CG(KSP ksp)542d71ae5a4SJacob Faibussowitsch PetscErrorCode KSPDestroy_CG(KSP ksp)
543d71ae5a4SJacob Faibussowitsch {
5444b9ad928SBarry Smith   KSP_CG *cg = (KSP_CG *)ksp->data;
5454b9ad928SBarry Smith 
5464b9ad928SBarry Smith   PetscFunctionBegin;
5479566063dSJacob Faibussowitsch   PetscCall(PetscFree4(cg->e, cg->d, cg->ee, cg->dd));
5489566063dSJacob Faibussowitsch   PetscCall(KSPDestroyDefault(ksp));
549fb01098fSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetObjectiveTarget_C", NULL));
5501a6f6928SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetRadius_C", NULL));
5519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetType_C", NULL));
5529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGUseSingleReduction_C", NULL));
553df8b9251SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGGetObjFcn_C", NULL));
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5554b9ad928SBarry Smith }
5564b9ad928SBarry Smith 
5574b9ad928SBarry Smith /*
55863017de6SPatrick Sanan      KSPView_CG - Prints information about the current Krylov method being used.
55963017de6SPatrick Sanan                   If your Krylov method has special options or flags that information
56063017de6SPatrick Sanan                   should be printed here.
5614b9ad928SBarry Smith */
KSPView_CG(KSP ksp,PetscViewer viewer)562d71ae5a4SJacob Faibussowitsch PetscErrorCode KSPView_CG(KSP ksp, PetscViewer viewer)
563d71ae5a4SJacob Faibussowitsch {
5644b9ad928SBarry Smith   KSP_CG   *cg = (KSP_CG *)ksp->data;
5659f196a02SMartin Diehl   PetscBool isascii;
5664b9ad928SBarry Smith 
5674b9ad928SBarry Smith   PetscFunctionBegin;
5689f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
5699f196a02SMartin Diehl   if (isascii) {
57030045a72SPatrick Sanan #if defined(PETSC_USE_COMPLEX)
5719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  variant %s\n", KSPCGTypes[cg->type]));
5724b9ad928SBarry Smith #endif
57348a46eb9SPierre Jolivet     if (cg->singlereduction) PetscCall(PetscViewerASCIIPrintf(viewer, "  using single-reduction variant\n"));
57430045a72SPatrick Sanan   }
5753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5764b9ad928SBarry Smith }
5774b9ad928SBarry Smith 
5784b9ad928SBarry Smith /*
5794b9ad928SBarry Smith     KSPSetFromOptions_CG - Checks the options database for options related to the
5804b9ad928SBarry Smith                            conjugate gradient method.
5814b9ad928SBarry Smith */
KSPSetFromOptions_CG(KSP ksp,PetscOptionItems PetscOptionsObject)582ce78bad3SBarry Smith PetscErrorCode KSPSetFromOptions_CG(KSP ksp, PetscOptionItems PetscOptionsObject)
583d71ae5a4SJacob Faibussowitsch {
584f0fe08ceSSatish Balay   KSP_CG   *cg = (KSP_CG *)ksp->data;
585c1bfdc7cSJed Brown   PetscBool flg;
5864b9ad928SBarry Smith 
5874b9ad928SBarry Smith   PetscFunctionBegin;
588d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "KSP CG and CGNE options");
589608df7e2SBarry Smith #if defined(PETSC_USE_COMPLEX)
5909371c9d4SSatish Balay   PetscCall(PetscOptionsEnum("-ksp_cg_type", "Matrix is Hermitian or complex symmetric", "KSPCGSetType", KSPCGTypes, (PetscEnum)cg->type, (PetscEnum *)&cg->type, NULL));
5914b9ad928SBarry Smith #endif
5929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ksp_cg_single_reduction", "Merge inner products into single MPI_Allreduce()", "KSPCGUseSingleReduction", cg->singlereduction, &cg->singlereduction, &flg));
5939566063dSJacob Faibussowitsch   if (flg) PetscCall(KSPCGUseSingleReduction(ksp, cg->singlereduction));
594d0609cedSBarry Smith   PetscOptionsHeadEnd();
5953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5964b9ad928SBarry Smith }
5974b9ad928SBarry Smith 
5984b9ad928SBarry Smith /*
5994b9ad928SBarry Smith     KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
6004b9ad928SBarry Smith                       This routine is registered below in KSPCreate_CG() and called from the
6014b9ad928SBarry Smith                       routine KSPCGSetType() (see the file cgtype.c).
6024b9ad928SBarry Smith */
KSPCGSetType_CG(KSP ksp,KSPCGType type)603d71ae5a4SJacob Faibussowitsch PetscErrorCode KSPCGSetType_CG(KSP ksp, KSPCGType type)
604d71ae5a4SJacob Faibussowitsch {
605608df7e2SBarry Smith   KSP_CG *cg = (KSP_CG *)ksp->data;
6064b9ad928SBarry Smith 
6074b9ad928SBarry Smith   PetscFunctionBegin;
6084b9ad928SBarry Smith   cg->type = type;
6093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6104b9ad928SBarry Smith }
6114b9ad928SBarry Smith 
61214579225SPatrick Sanan /*
61314579225SPatrick Sanan     KSPCGUseSingleReduction_CG
61414579225SPatrick Sanan 
61514579225SPatrick Sanan     This routine sets a flag to use a variant of CG. Note that (in somewhat
61614579225SPatrick Sanan     atypical fashion) it also swaps out the routine called when KSPSolve()
61714579225SPatrick Sanan     is invoked.
61814579225SPatrick Sanan */
KSPCGUseSingleReduction_CG(KSP ksp,PetscBool flg)619d71ae5a4SJacob Faibussowitsch static PetscErrorCode KSPCGUseSingleReduction_CG(KSP ksp, PetscBool flg)
620d71ae5a4SJacob Faibussowitsch {
621608df7e2SBarry Smith   KSP_CG *cg = (KSP_CG *)ksp->data;
622608df7e2SBarry Smith 
623608df7e2SBarry Smith   PetscFunctionBegin;
624608df7e2SBarry Smith   cg->singlereduction = flg;
62514579225SPatrick Sanan   if (cg->singlereduction) {
6268d36b627SPatrick Sanan     ksp->ops->solve = KSPSolve_CG_SingleReduction;
627df8b9251SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGGetObjFcn_C", NULL));
62814579225SPatrick Sanan   } else {
62914579225SPatrick Sanan     ksp->ops->solve = KSPSolve_CG;
630df8b9251SStefano Zampini     PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGGetObjFcn_C", KSPCGGetObjFcn_CG));
63114579225SPatrick Sanan   }
6323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
633608df7e2SBarry Smith }
634608df7e2SBarry Smith 
KSPBuildResidual_CG(KSP ksp,Vec t,Vec v,Vec * V)635d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode KSPBuildResidual_CG(KSP ksp, Vec t, Vec v, Vec *V)
636d71ae5a4SJacob Faibussowitsch {
6378258c059SStefano Zampini   PetscFunctionBegin;
6389566063dSJacob Faibussowitsch   PetscCall(VecCopy(ksp->work[0], v));
6398258c059SStefano Zampini   *V = v;
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6418258c059SStefano Zampini }
6428258c059SStefano Zampini 
6438742991fSBarry Smith /*MC
6440b4b7b1cSBarry Smith    KSPCG - The Preconditioned Conjugate Gradient (PCG) iterative method {cite}`hs:52` and {cite}`malek2014preconditioning` for solving linear systems using `KSP`.
6458742991fSBarry Smith 
6468742991fSBarry Smith    Options Database Keys:
647dc80d9caSBarry Smith +   -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian, see `KSPCGSetType()`
64806e420deSBarry Smith .   -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric
649dc80d9caSBarry Smith -   -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single `MPI_Allreduce()` call, see `KSPCGUseSingleReduction()`
6508742991fSBarry Smith 
6518742991fSBarry Smith    Level: beginner
6528742991fSBarry Smith 
65395452b02SPatrick Sanan    Notes:
6540b4b7b1cSBarry Smith    The `KSPCG` method requires both the matrix and preconditioner to be symmetric positive (or negative) (semi) definite.
6558742991fSBarry Smith 
6560b4b7b1cSBarry Smith    `KSPCG` is the best Krylov method, `KSPType`, when the matrix and preconditioner are symmetric positive definite (SPD).
6570b4b7b1cSBarry Smith 
6580b4b7b1cSBarry Smith    Only left preconditioning is supported with `KSPCG`; there are several ways to motivate preconditioned CG, but they all produce the same algorithm.
6590b4b7b1cSBarry Smith    One can interpret preconditioning $A$ with $B$ to mean any of the following\:
660562efe2eSBarry Smith .vb
6610b4b7b1cSBarry Smith    (1) Solve a left-preconditioned system $BAx = Bb $, using $ B^{-1}$ to define an inner product in the algorithm.
6620b4b7b1cSBarry Smith    (2) Solve a right-preconditioned system $ABy = b, x = By,$ using $B$ to define an inner product in the algorithm.
6630b4b7b1cSBarry Smith    (3) Solve a symmetrically-preconditioned system, $ E^TAEy = E^Tb, x = Ey, $ where $B = EE^T.$
6640b4b7b1cSBarry Smith    (4) Solve $Ax=b$ with CG, but use the inner product defined by $B$ to define the method.
665789736e1SBarry Smith    In all cases, the resulting algorithm only requires application of $B$ to vectors, the other inner-product does not appear explicitly in the code
666562efe2eSBarry Smith .ve
66737c1d77aSPatrick Sanan 
66837c1d77aSPatrick Sanan    For complex numbers there are two different CG methods, one for Hermitian symmetric matrices and one for non-Hermitian symmetric matrices. Use
669dc80d9caSBarry Smith    `KSPCGSetType()` to indicate which type you are using.
670dc80d9caSBarry Smith 
671dc80d9caSBarry Smith    One can use `KSPSetComputeEigenvalues()` and `KSPComputeEigenvalues()` to compute the eigenvalues of the (preconditioned) operator
672828413b8SBarry Smith 
6730b4b7b1cSBarry Smith    There are two pipelined implementations of CG in PETSc `KSPPIPECG` and `KSPGROPPCG`. These may perform better for very large
6740b4b7b1cSBarry Smith    numbers of MPI processes since they overlap communication and computation so the reduction operations in CG, that is inner products and norms,
6750b4b7b1cSBarry Smith    do not dominate the compute time.
6760b4b7b1cSBarry Smith 
677562efe2eSBarry Smith    Developer Note:
6780b4b7b1cSBarry Smith    KSPSolve_CG() should actually query the matrix to determine if it is Hermitian or symmetric and NOT require the user to
679dc80d9caSBarry Smith    indicate it to the `KSP` object.
680828413b8SBarry Smith 
6811cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPSetComputeEigenvalues()`, `KSPComputeEigenvalues()`
682db781477SPatrick Sanan           `KSPCGSetType()`, `KSPCGUseSingleReduction()`, `KSPPIPECG`, `KSPGROPPCG`
6838742991fSBarry Smith M*/
684dc80d9caSBarry Smith 
685dc80d9caSBarry Smith /*
686dc80d9caSBarry Smith     KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
687dc80d9caSBarry Smith        function pointers for all the routines it needs to call (KSPSolve_CG() etc)
688dc80d9caSBarry Smith 
689dc80d9caSBarry Smith     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
690dc80d9caSBarry Smith */
KSPCreate_CG(KSP ksp)691d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode KSPCreate_CG(KSP ksp)
692d71ae5a4SJacob Faibussowitsch {
6934b9ad928SBarry Smith   KSP_CG *cg;
6944b9ad928SBarry Smith 
6954b9ad928SBarry Smith   PetscFunctionBegin;
6964dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&cg));
697*fc2fb351SPierre Jolivet   cg->type    = !PetscDefined(USE_COMPLEX) ? KSP_CG_SYMMETRIC : KSP_CG_HERMITIAN;
698fb01098fSStefano Zampini   cg->obj_min = 0.0;
6994b9ad928SBarry Smith   ksp->data   = (void *)cg;
7009e568555SJed Brown 
7019566063dSJacob Faibussowitsch   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 3));
7029566063dSJacob Faibussowitsch   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2));
7039566063dSJacob Faibussowitsch   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2));
7049566063dSJacob Faibussowitsch   PetscCall(KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 1));
7054b9ad928SBarry Smith 
7064b9ad928SBarry Smith   /*
7074b9ad928SBarry Smith        Sets the functions that are associated with this data structure
7084b9ad928SBarry Smith        (in C++ this is the same as defining virtual functions)
7094b9ad928SBarry Smith   */
7104b9ad928SBarry Smith   ksp->ops->setup          = KSPSetUp_CG;
7114b9ad928SBarry Smith   ksp->ops->solve          = KSPSolve_CG;
7124b9ad928SBarry Smith   ksp->ops->destroy        = KSPDestroy_CG;
7134b9ad928SBarry Smith   ksp->ops->view           = KSPView_CG;
7144b9ad928SBarry Smith   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
715fa0ddf94SBarry Smith   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
7168258c059SStefano Zampini   ksp->ops->buildresidual  = KSPBuildResidual_CG;
7174b9ad928SBarry Smith 
7184b9ad928SBarry Smith   /*
7194b9ad928SBarry Smith       Attach the function KSPCGSetType_CG() to this object. The routine
7204b9ad928SBarry Smith       KSPCGSetType() checks for this attached function and calls it if it finds
7214b9ad928SBarry Smith       it. (Sort of like a dynamic member function that can be added at run time
7224b9ad928SBarry Smith   */
7239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetType_C", KSPCGSetType_CG));
7249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGUseSingleReduction_C", KSPCGUseSingleReduction_CG));
7251a6f6928SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetRadius_C", KSPCGSetRadius_CG));
726fb01098fSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGSetObjectiveTarget_C", KSPCGSetObjectiveTarget_CG));
727df8b9251SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)ksp, "KSPCGGetObjFcn_C", KSPCGGetObjFcn_CG));
7283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7294b9ad928SBarry Smith }
730