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: cgimpl.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/cgimpl.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,nwork = 3;
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: if (cgP->singlereduction) nwork += 2;
76: KSPDefaultGetWork(ksp,nwork);
78: /*
79: If user requested computations of eigenvalues then allocate work
80: work space needed
81: */
82: if (ksp->calc_sings) {
83: /* get space to store tridiagonal matrix for Lanczos */
84: PetscMalloc4(maxit+1,PetscScalar,&cgP->e,maxit+1,PetscScalar,&cgP->d,maxit+1,PetscReal,&cgP->ee,maxit+1,PetscReal,&cgP->dd);
85: PetscLogObjectMemory(ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));
86: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
87: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_CG;
88: }
89: return(0);
90: }
92: /*
93: KSPSolve_CG - This routine actually applies the conjugate gradient method
95: This routine is MUCH too messy. I has too many options (norm type and single reduction) embedded making the code confusing and likely to be buggy.
97: Input Parameter:
98: . ksp - the Krylov space object that was set to use conjugate gradient, by, for
99: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
100: */
103: PetscErrorCode KSPSolve_CG(KSP ksp)
104: {
106: PetscInt i,stored_max_it,eigs;
107: PetscScalar dpi = 0.0,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0,delta,dpiold;
108: PetscReal dp = 0.0;
109: Vec X,B,Z,R,P,S,W;
110: KSP_CG *cg;
111: Mat Amat,Pmat;
112: MatStructure pflag;
113: PetscTruth diagonalscale;
116: PCDiagonalScale(ksp->pc,&diagonalscale);
117: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
119: cg = (KSP_CG*)ksp->data;
120: eigs = ksp->calc_sings;
121: stored_max_it = ksp->max_it;
122: X = ksp->vec_sol;
123: B = ksp->vec_rhs;
124: R = ksp->work[0];
125: Z = ksp->work[1];
126: P = ksp->work[2];
127: if (cg->singlereduction) {
128: S = ksp->work[3];
129: W = ksp->work[4];
130: } else {
131: S = 0; /* unused */
132: W = Z;
133: }
135: #if !defined(PETSC_USE_COMPLEX)
136: #define VecXDot(x,y,a) VecDot(x,y,a)
137: #else
138: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
139: #endif
141: if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
142: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
144: ksp->its = 0;
145: if (!ksp->guess_zero) {
146: KSP_MatMult(ksp,Amat,X,R); /* r <- b - Ax */
147: VecAYPX(R,-1.0,B);
148: } else {
149: VecCopy(B,R); /* r <- b (x is 0) */
150: }
152: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
153: KSP_PCApply(ksp,R,Z); /* z <- Br */
154: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z = e'*A'*B'*B*A'*e' */
155: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
156: VecNorm(R,NORM_2,&dp); /* dp <- r'*r = e'*A'*A*e */
157: } else if (ksp->normtype == KSP_NORM_NATURAL) {
158: KSP_PCApply(ksp,R,Z); /* z <- Br */
159: if (cg->singlereduction) {
160: KSP_MatMult(ksp,Amat,Z,S);
161: VecXDot(Z,S,&delta);
162: }
163: VecXDot(Z,R,&beta); /* beta <- z'*r */
164: if PetscIsInfOrNanScalar(beta) SETERRQ(PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
165: dp = sqrt(PetscAbsScalar(beta)); /* dp <- r'*z = r'*B*r = e'*A'*B*A*e */
166: } else dp = 0.0;
167: KSPLogResidualHistory(ksp,dp);
168: KSPMonitor(ksp,0,dp); /* call any registered monitor routines */
169: ksp->rnorm = dp;
171: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
172: if (ksp->reason) return(0);
174: if (ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)){
175: KSP_PCApply(ksp,R,Z); /* z <- Br */
176: }
177: if (ksp->normtype != KSP_NORM_NATURAL){
178: if (cg->singlereduction) {
179: KSP_MatMult(ksp,Amat,Z,S);
180: VecXDot(Z,S,&delta);
181: }
182: VecXDot(Z,R,&beta); /* beta <- z'*r */
183: if PetscIsInfOrNanScalar(beta) SETERRQ(PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
184: }
185:
186: i = 0;
187: do {
188: ksp->its = i+1;
189: if (beta == 0.0) {
190: ksp->reason = KSP_CONVERGED_ATOL;
191: PetscInfo(ksp,"converged due to beta = 0\n");
192: break;
193: #if !defined(PETSC_USE_COMPLEX)
194: } else if (beta < 0.0) {
195: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
196: PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
197: break;
198: #endif
199: }
200: if (!i) {
201: VecCopy(Z,P); /* p <- z */
202: b = 0.0;
203: } else {
204: b = beta/betaold;
205: if (eigs) {
206: if (ksp->max_it != stored_max_it) {
207: SETERRQ(PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
208: }
209: e[i] = sqrt(PetscAbsScalar(b))/a;
210: }
211: VecAYPX(P,b,Z); /* p <- z + b* p */
212: }
213: dpiold = dpi;
214: if (!cg->singlereduction || !i) {
215: KSP_MatMult(ksp,Amat,P,W); /* w <- Kp */
216: VecXDot(P,W,&dpi); /* dpi <- p'w */
217: } else {
218: VecAYPX(W,beta/betaold,S); /* w <- Kp */
219: dpi = delta - beta*beta*dpiold/(betaold*betaold); /* dpi <- p'w */
220: }
221: betaold = beta;
222: if PetscIsInfOrNanScalar(dpi) SETERRQ(PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
224: if (PetscRealPart(dpi) <= 0.0) {
225: ksp->reason = KSP_DIVERGED_INDEFINITE_MAT;
226: PetscInfo(ksp,"diverging due to indefinite or negative definite matrix\n");
227: break;
228: }
229: a = beta/dpi; /* a = beta/p'w */
230: if (eigs) {
231: d[i] = sqrt(PetscAbsScalar(b))*e[i] + 1.0/a;
232: }
233: VecAXPY(X,a,P); /* x <- x + ap */
234: VecAXPY(R,-a,W); /* r <- r - aw */
235: if (ksp->normtype == KSP_NORM_PRECONDITIONED && ksp->chknorm < i+2) {
236: KSP_PCApply(ksp,R,Z); /* z <- Br */
237: if (cg->singlereduction) {
238: KSP_MatMult(ksp,Amat,Z,S);
239: }
240: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
241: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED && ksp->chknorm < i+2) {
242: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
243: } else if (ksp->normtype == KSP_NORM_NATURAL) {
244: KSP_PCApply(ksp,R,Z); /* z <- Br */
245: if (cg->singlereduction) {
246: PetscScalar tmp[2];
247: Vec vecs[2];
248: vecs[0] = S; vecs[1] = R;
249: KSP_MatMult(ksp,Amat,Z,S);
250: /*VecXDot(Z,S,&delta);
251: VecXDot(Z,R,&beta); */ /* beta <- r'*z */
252: VecMDot(Z,2,vecs,tmp);
253: delta = tmp[0]; beta = tmp[1];
254: } else {
255: VecXDot(Z,R,&beta); /* beta <- r'*z */
256: }
257: if PetscIsInfOrNanScalar(beta) SETERRQ(PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
258: dp = sqrt(PetscAbsScalar(beta));
259: } else {
260: dp = 0.0;
261: }
262: ksp->rnorm = dp;
263: KSPLogResidualHistory(ksp,dp);
264: KSPMonitor(ksp,i+1,dp);
265: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
266: if (ksp->reason) break;
268: if ((ksp->normtype != KSP_NORM_PRECONDITIONED && (ksp->normtype != KSP_NORM_NATURAL)) || (ksp->chknorm >= i+2)){
269: KSP_PCApply(ksp,R,Z); /* z <- Br */
270: if (cg->singlereduction) {
271: KSP_MatMult(ksp,Amat,Z,S);
272: }
273: }
274: if ((ksp->normtype != KSP_NORM_NATURAL) || (ksp->chknorm >= i+2)){
275: if (cg->singlereduction) {
276: PetscScalar tmp[2];
277: Vec vecs[2];
278: vecs[0] = S; vecs[1] = R;
279: /* VecXDot(Z,R,&beta); */ /* beta <- z'*r */
280: /* VecXDot(Z,S,&delta);*/
281: VecMDot(Z,2,vecs,tmp);
282: delta = tmp[0]; beta = tmp[1];
283: } else {
284: VecXDot(Z,R,&beta); /* beta <- z'*r */
285: }
286: if PetscIsInfOrNanScalar(beta) SETERRQ(PETSC_ERR_FP,"Infinite or not-a-number generated in dot product");
287: }
289: i++;
290: } while (i<ksp->max_it);
291: if (i >= ksp->max_it) {
292: ksp->reason = KSP_DIVERGED_ITS;
293: }
294: return(0);
295: }
296: /*
297: KSPDestroy_CG - Frees all memory space used by the Krylov method
299: */
302: PetscErrorCode KSPDestroy_CG(KSP ksp)
303: {
304: KSP_CG *cg = (KSP_CG*)ksp->data;
308: /* free space used for singular value calculations */
309: if (ksp->calc_sings) {
310: PetscFree4(cg->e,cg->dd,cg->ee,cg->dd);
311: }
312: KSPDefaultDestroy(ksp);
313: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","",PETSC_NULL);
314: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGUseSingleReduction_C","",PETSC_NULL);
315: return(0);
316: }
318: /*
319: KSPView_CG - Prints information about the current Krylov method being used
321: Currently this only prints information to a file (or stdout) about the
322: symmetry of the problem. If your Krylov method has special options or
323: flags that information should be printed here.
325: */
328: PetscErrorCode KSPView_CG(KSP ksp,PetscViewer viewer)
329: {
330: #if defined(PETSC_USE_COMPLEX)
331: KSP_CG *cg = (KSP_CG *)ksp->data;
333: PetscTruth iascii;
336: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
337: if (iascii) {
338: PetscViewerASCIIPrintf(viewer," CG or CGNE: variant %s\n",KSPCGTypes[cg->type]);
339: } else {
340: SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
341: }
342: #endif
343: return(0);
344: }
346: /*
347: KSPSetFromOptions_CG - Checks the options database for options related to the
348: conjugate gradient method.
349: */
352: PetscErrorCode KSPSetFromOptions_CG(KSP ksp)
353: {
355: KSP_CG *cg = (KSP_CG *)ksp->data;
358: PetscOptionsHead("KSP CG and CGNE options");
359: #if defined(PETSC_USE_COMPLEX)
360: PetscOptionsEnum("-ksp_cg_type","Matrix is Hermitian or complex symmetric","KSPCGSetType",KSPCGTypes,(PetscEnum)cg->type,
361: (PetscEnum*)&cg->type,PETSC_NULL);
362: #endif
363: PetscOptionsTruth("-ksp_cg_single_reduction","Merge inner products into single MPI_Allreduce()",
364: "KSPCGUseSingleReduction",cg->singlereduction,&cg->singlereduction,PETSC_NULL);
365: PetscOptionsTail();
366: return(0);
367: }
369: /*
370: KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
371: This routine is registered below in KSPCreate_CG() and called from the
372: routine KSPCGSetType() (see the file cgtype.c).
375: */
379: PetscErrorCode KSPCGSetType_CG(KSP ksp,KSPCGType type)
380: {
381: KSP_CG *cg = (KSP_CG *)ksp->data;
384: cg->type = type;
385: return(0);
386: }
392: PetscErrorCode KSPCGUseSingleReduction_CG(KSP ksp,PetscTruth flg)
393: {
394: KSP_CG *cg = (KSP_CG *)ksp->data;
397: cg->singlereduction = flg;
398: return(0);
399: }
402: /*
403: KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
404: function pointers for all the routines it needs to call (KSPSolve_CG() etc)
407: */
408: /*MC
409: KSPCG - The preconditioned conjugate gradient (PCG) iterative method
411: Options Database Keys:
412: + -ksp_cg_type Hermitian - (for complex matrices only) indicates the matrix is Hermitian
413: . -ksp_cg_type symmetric - (for complex matrices only) indicates the matrix is symmetric
414: - -ksp_cg_single_reduction - performs both inner products needed in the algorithm with a single MPI_Allreduce() call, see KSPCGUseSingleReduction()
416: Level: beginner
418: Notes: The PCG method requires both the matrix and preconditioner to
419: be symmetric positive (semi) definite
421: References:
422: Methods of Conjugate Gradients for Solving Linear Systems, Magnus R. Hestenes and Eduard Stiefel,
423: Journal of Research of the National Bureau of Standards Vol. 49, No. 6, December 1952 Research Paper 2379
424: pp. 409--436.
426: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
427: KSPCGSetType(), KSPCGUseSingleReduction()
429: M*/
433: PetscErrorCode KSPCreate_CG(KSP ksp)
434: {
436: KSP_CG *cg;
439: PetscNewLog(ksp,KSP_CG,&cg);
440: #if !defined(PETSC_USE_COMPLEX)
441: cg->type = KSP_CG_SYMMETRIC;
442: #else
443: cg->type = KSP_CG_HERMITIAN;
444: #endif
445: ksp->data = (void*)cg;
446: ksp->pc_side = PC_LEFT;
448: /*
449: Sets the functions that are associated with this data structure
450: (in C++ this is the same as defining virtual functions)
451: */
452: ksp->ops->setup = KSPSetUp_CG;
453: ksp->ops->solve = KSPSolve_CG;
454: ksp->ops->destroy = KSPDestroy_CG;
455: ksp->ops->view = KSPView_CG;
456: ksp->ops->setfromoptions = KSPSetFromOptions_CG;
457: ksp->ops->buildsolution = KSPDefaultBuildSolution;
458: ksp->ops->buildresidual = KSPDefaultBuildResidual;
460: /*
461: Attach the function KSPCGSetType_CG() to this object. The routine
462: KSPCGSetType() checks for this attached function and calls it if it finds
463: it. (Sort of like a dynamic member function that can be added at run time
464: */
465: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG", KSPCGSetType_CG);
466: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGUseSingleReduction_C","KSPCGUseSingleReduction_CG", KSPCGUseSingleReduction_CG);
467: return(0);
468: }