Actual source code: cg.c
1: #define PETSCKSP_DLL
3: /*
4: This file implements the conjugate gradient method in PETSc as part of
5: KSP. You can use this as a starting point for implementing your own
6: Krylov method that is not provided with PETSc.
8: The following basic routines are required for each Krylov method.
9: KSPCreate_XXX() - Creates the Krylov context
10: KSPSetFromOptions_XXX() - Sets runtime options
11: KSPSolve_XXX() - Runs the Krylov method
12: KSPDestroy_XXX() - Destroys the Krylov context, freeing all
13: memory it needed
14: Here the "_XXX" denotes a particular implementation, in this case
15: we use _CG (e.g. KSPCreate_CG, KSPDestroy_CG). These routines are
16: are actually called vai the common user interface routines
17: KSPSetType(), KSPSetFromOptions(), KSPSolve(), and KSPDestroy() so the
18: application code interface remains identical for all preconditioners.
20: Other basic routines for the KSP objects include
21: KSPSetUp_XXX()
22: KSPView_XXX() - Prints details of solver being used.
24: Detailed notes:
25: By default, this code implements the CG (Conjugate Gradient) method,
26: which is valid for real symmetric (and complex Hermitian) positive
27: definite matrices. Note that for the complex Hermitian case, the
28: VecDot() arguments within the code MUST remain in the order given
29: for correct computation of inner products.
31: Reference: Hestenes and Steifel, 1952.
33: By switching to the indefinite vector inner product, VecTDot(), the
34: same code is used for the complex symmetric case as well. The user
35: must call KSPCGSetType(ksp,KSP_CG_SYMMETRIC) or use the option
36: -ksp_cg_type symmetric to invoke this variant for the complex case.
37: Note, however, that the complex symmetric code is NOT valid for
38: all such matrices ... and thus we don't recommend using this method.
39: */
40: /*
41: cgctx.h defines the simple data structured used to store information
42: related to the type of matrix (e.g. complex symmetric) being solved and
43: data used during the optional Lanczo process used to compute eigenvalues
44: */
45: #include src/ksp/ksp/impls/cg/cgctx.h
46: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
47: EXTERN PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);
49: /*
50: KSPSetUp_CG - Sets up the workspace needed by the CG method.
52: This is called once, usually automatically by KSPSolve() or KSPSetUp()
53: but can be called directly by KSPSetUp()
54: */
57: PetscErrorCode KSPSetUp_CG(KSP ksp)
58: {
59: KSP_CG *cgP = (KSP_CG*)ksp->data;
61: PetscInt maxit = ksp->max_it;
64: /*
65: This implementation of CG only handles left preconditioning
66: so generate an error otherwise.
67: */
68: if (ksp->pc_side == PC_RIGHT) {
69: SETERRQ(PETSC_ERR_SUP,"No right preconditioning for KSPCG");
70: } else if (ksp->pc_side == PC_SYMMETRIC) {
71: SETERRQ(PETSC_ERR_SUP,"No symmetric preconditioning for KSPCG");
72: }
74: /* get work vectors needed by CG */
75: KSPDefaultGetWork(ksp,3);
77: /*
78: If user requested computations of eigenvalues then allocate work
79: work space needed
80: */
81: if (ksp->calc_sings) {
82: /* get space to store tridiagonal matrix for Lanczos */
83: PetscMalloc(2*(maxit+1)*sizeof(PetscScalar),&cgP->e);
84: PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
85: cgP->d = cgP->e + maxit + 1;
86: PetscMalloc(2*(maxit+1)*sizeof(PetscReal),&cgP->ee);
87: PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
88: cgP->dd = cgP->ee + maxit + 1;
89: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
90: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_CG;
91: }
92: return(0);
93: }
95: /*
96: KSPSolve_CG - This routine actually applies the conjugate gradient
97: method
99: Input Parameter:
100: . ksp - the Krylov space object that was set to use conjugate gradient, by, for
101: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
102: */
105: PetscErrorCode KSPSolve_CG(KSP ksp)
106: {
108: PetscInt i,stored_max_it,eigs;
109: PetscScalar dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0;
110: PetscReal dp = 0.0;
111: Vec X,B,Z,R,P;
112: KSP_CG *cg;
113: Mat Amat,Pmat;
114: MatStructure pflag;
115: PetscTruth diagonalscale;
118: PCDiagonalScale(ksp->pc,&diagonalscale);
119: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
121: cg = (KSP_CG*)ksp->data;
122: eigs = ksp->calc_sings;
123: stored_max_it = ksp->max_it;
124: X = ksp->vec_sol;
125: B = ksp->vec_rhs;
126: R = ksp->work[0];
127: Z = ksp->work[1];
128: P = ksp->work[2];
130: #if !defined(PETSC_USE_COMPLEX)
131: #define VecXDot(x,y,a) VecDot(x,y,a)
132: #else
133: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
134: #endif
136: if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
137: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
139: ksp->its = 0;
140: if (!ksp->guess_zero) {
141: KSP_MatMult(ksp,Amat,X,R); /* r <- b - Ax */
142: VecAYPX(R,-1.0,B);
143: } else {
144: VecCopy(B,R); /* r <- b (x is 0) */
145: }
147: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
148: KSP_PCApply(ksp,R,Z); /* z <- Br */
149: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z = e'*A'*B'*B*A'*e' */
150: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
151: VecNorm(R,NORM_2,&dp); /* dp <- r'*r = e'*A'*A*e */
152: } else if (ksp->normtype == KSP_NORM_NATURAL) {
153: KSP_PCApply(ksp,R,Z); /* z <- Br */
154: VecXDot(Z,R,&beta);CHKFPQ(beta); /* beta <- z'*r */
155: dp = sqrt(PetscAbsScalar(beta)); /* dp <- r'*z = r'*B*r = e'*A'*B*A*e */
156: } else dp = 0.0;
157: KSPLogResidualHistory(ksp,dp);
158: KSPMonitor(ksp,0,dp); /* call any registered monitor routines */
159: ksp->rnorm = dp;
161: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
162: if (ksp->reason) return(0);
164: if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)){
165: KSP_PCApply(ksp,R,Z); /* z <- Br */
166: }
167: if (ksp->normtype != KSP_NORM_NATURAL){
168: VecXDot(Z,R,&beta);CHKFPQ(beta); /* beta <- z'*r */
169: }
171: i = 0;
172: do {
173: ksp->its = i+1;
174: if (beta == 0.0) {
175: ksp->reason = KSP_CONVERGED_ATOL;
176: PetscInfo(ksp,"converged due to beta = 0\n");
177: break;
178: #if !defined(PETSC_USE_COMPLEX)
179: } else if (beta < 0.0) {
180: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
181: PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
182: break;
183: #endif
184: }
185: if (!i) {
186: VecCopy(Z,P); /* p <- z */
187: b = 0.0;
188: } else {
189: b = beta/betaold;
190: if (eigs) {
191: if (ksp->max_it != stored_max_it) {
192: SETERRQ(PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
193: }
194: e[i] = sqrt(PetscAbsScalar(b))/a;
195: }
196: VecAYPX(P,b,Z); /* p <- z + b* p */
197: }
198: betaold = beta;
199: KSP_MatMult(ksp,Amat,P,Z); /* z <- Kp */
200: VecXDot(P,Z,&dpi);CHKFPQ(dpi); /* dpi <- z'p */
202: if (PetscRealPart(dpi) <= 0.0) {
203: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
204: PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
205: break;
206: }
207: a = beta/dpi; /* a = beta/p'z */
208: if (eigs) {
209: d[i] = sqrt(PetscAbsScalar(b))*e[i] + 1.0/a;
210: }
211: VecAXPY(X,a,P); /* x <- x + ap */
212: VecAXPY(R,-a,Z); /* r <- r - az */
213: if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
214: KSP_PCApply(ksp,R,Z); /* z <- Br */
215: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
216: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
217: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
218: } else if (ksp->normtype == KSP_NORM_NATURAL) {
219: KSP_PCApply(ksp,R,Z); /* z <- Br */
220: VecXDot(Z,R,&beta);CHKFPQ(beta); /* beta <- r'*z */
221: dp = sqrt(PetscAbsScalar(beta));
222: } else {
223: dp = 0.0;
224: }
225: ksp->rnorm = dp;
226: KSPLogResidualHistory(ksp,dp);
227: KSPMonitor(ksp,i+1,dp);
228: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
229: if (ksp->reason) break;
231: if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)){
232: KSP_PCApply(ksp,R,Z); /* z <- Br */
233: }
234: if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)){
235: VecXDot(Z,R,&beta);CHKFPQ(beta); /* beta <- z'*r */
236: }
238: i++;
239: } while (i<ksp->max_it);
240: if (i >= ksp->max_it) {
241: ksp->reason = KSP_DIVERGED_ITS;
242: }
243: return(0);
244: }
245: /*
246: KSPDestroy_CG - Frees all memory space used by the Krylov method
248: */
251: PetscErrorCode KSPDestroy_CG(KSP ksp)
252: {
253: KSP_CG *cg = (KSP_CG*)ksp->data;
257: /* free space used for singular value calculations */
258: if (ksp->calc_sings) {
259: PetscFree(cg->e);
260: PetscFree(cg->ee);
261: }
262: KSPDefaultDestroy(ksp);
263: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","",PETSC_NULL);
264: return(0);
265: }
267: /*
268: KSPView_CG - Prints information about the current Krylov method being used
270: Currently this only prints information to a file (or stdout) about the
271: symmetry of the problem. If your Krylov method has special options or
272: flags that information should be printed here.
274: */
277: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
278: {
279: #if defined(PETSC_USE_COMPLEX)
280: KSP_CG *cg = (KSP_CG *)ksp->data;
282: PetscTruth iascii;
285: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
286: if (iascii) {
287: PetscViewerASCIIPrintf(viewer," CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
288: } else {
289: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
290: }
291: #endif
292: return(0);
293: }
295: /*
296: KSPSetFromOptions_CG - Checks the options database for options related to the
297: conjugate gradient method.
298: */
301: PetscErrorCode KSPSetFromOptions_CG(KSP ksp)
302: {
303: #if defined(PETSC_USE_COMPLEX)
305: KSP_CG *cg = (KSP_CG *)ksp->data;
306: #endif
309: #if defined(PETSC_USE_COMPLEX)
310: PetscOptionsHead("KSP CG and CGNE options");
311: PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
312: (PetscEnum*)&cg->type,PETSC_NULL);
313: PetscOptionsTail();
314: #endif
315: return(0);
316: }
318: /*
319: KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
320: This routine is registered below in KSPCreate_CG() and called from the
321: routine KSPCGSetType() (see the file cgtype.c).
324: */
328: PetscErrorCode KSPCGSetType_CG(KSP ksp,KSPCGType type)
329: {
330: KSP_CG *cg;
333: cg = (KSP_CG *)ksp->data;
334: cg->type = type;
335: return(0);
336: }
339: /*
340: KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
341: function pointers for all the routines it needs to call (KSPSolve_CG() etc)
344: */
345: /*MC
346: KSPCG - The preconditioned conjugate gradient (PCG) iterative method
348: Options Database Keys:
349: + -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian
350: - -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric
352: Level: beginner
354: Notes: The PCG method requires both the matrix and preconditioner to
355: be symmetric positive (semi) definite
357: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
358: KSPCGSetType()
360: M*/
364: PetscErrorCode KSPCreate_CG(KSP ksp)
365: {
367: KSP_CG *cg;
370: PetscNewLog(ksp,KSP_CG,&cg);
371: #if !defined(PETSC_USE_COMPLEX)
372: cg->type = KSP_CG_SYMMETRIC;
373: #else
374: cg->type = KSP_CG_HERMITIAN;
375: #endif
376: ksp->data = (void*)cg;
377: ksp->pc_side = PC_LEFT;
379: /*
380: Sets the functions that are associated with this data structure
381: (in C++ this is the same as defining virtual functions)
382: */
383: ksp->ops->setup = KSPSetUp_CG;
384: ksp->ops->solve = KSPSolve_CG;
385: ksp->ops->destroy = KSPDestroy_CG;
386: ksp->ops->view = KSPView_CG;
387: ksp->ops->setfromoptions = KSPSetFromOptions_CG;
388: ksp->ops->buildsolution = KSPDefaultBuildSolution;
389: ksp->ops->buildresidual = KSPDefaultBuildResidual;
391: /*
392: Attach the function KSPCGSetType_CG() to this object. The routine
393: KSPCGSetType() checks for this attached function and calls it if it finds
394: it. (Sort of like a dynamic member function that can be added at run time
395: */
396: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C",
397: "KSPCGSetType_CG",
398: KSPCGSetType_CG);
399: return(0);
400: }