Actual source code: cg.c
1: /*$Id: cg.c,v 1.117 2001/08/07 03:03:50 balay Exp $*/
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_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 int KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
47: EXTERN int KSPComputeEigenvalues_CG(KSP,int,PetscReal *,PetscReal *,int *);
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: int KSPSetUp_CG(KSP ksp)
58: {
59: KSP_CG *cgP = (KSP_CG*)ksp->data;
60: int maxit = ksp->max_it,ierr;
63: /*
64: This implementation of CG only handles left preconditioning
65: so generate an error otherwise.
66: */
67: if (ksp->pc_side == PC_RIGHT) {
68: SETERRQ(2,"No right preconditioning for KSPCG");
69: } else if (ksp->pc_side == PC_SYMMETRIC) {
70: SETERRQ(2,"No symmetric preconditioning for KSPCG");
71: }
73: /* get work vectors needed by CG */
74: KSPDefaultGetWork(ksp,3);
76: /*
77: If user requested computations of eigenvalues then allocate work
78: work space needed
79: */
80: if (ksp->calc_sings) {
81: /* get space to store tridiagonal matrix for Lanczos */
82: PetscMalloc(2*(maxit+1)*sizeof(PetscScalar),&cgP->e);
83: PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
84: cgP->d = cgP->e + maxit + 1;
85: PetscMalloc(2*(maxit+1)*sizeof(PetscReal),&cgP->ee);
86: PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
87: cgP->dd = cgP->ee + maxit + 1;
88: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
89: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_CG;
90: }
91: return(0);
92: }
94: /*
95: KSPSolve_CG - This routine actually applies the conjugate gradient
96: method
98: Input Parameter:
99: . ksp - the Krylov space object that was set to use conjugate gradient, by, for
100: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
102: Output Parameter:
103: . its - number of iterations used
105: */
108: int KSPSolve_CG(KSP ksp)
109: {
110: int ierr,i,stored_max_it,eigs;
111: PetscScalar dpi,a = 1.0,beta,betaold = 1.0,b,*e = 0,*d = 0,mone = -1.0,ma;
112: PetscReal dp = 0.0;
113: Vec X,B,Z,R,P;
114: KSP_CG *cg;
115: Mat Amat,Pmat;
116: MatStructure pflag;
117: PetscTruth diagonalscale;
120: PCDiagonalScale(ksp->B,&diagonalscale);
121: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
123: cg = (KSP_CG*)ksp->data;
124: eigs = ksp->calc_sings;
125: stored_max_it = ksp->max_it;
126: X = ksp->vec_sol;
127: B = ksp->vec_rhs;
128: R = ksp->work[0];
129: Z = ksp->work[1];
130: P = ksp->work[2];
132: #if !defined(PETSC_USE_COMPLEX)
133: #define VecXDot(x,y,a) VecDot(x,y,a)
134: #else
135: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
136: #endif
138: if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; b = 0.0; }
139: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
141: ksp->its = 0;
142: if (!ksp->guess_zero) {
143: KSP_MatMult(ksp,Amat,X,R); /* r <- b - Ax */
144: VecAYPX(&mone,B,R);
145: } else {
146: VecCopy(B,R); /* r <- b (x is 0) */
147: }
148: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- Br */
149: VecXDot(Z,R,&beta);
150: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
151: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
152: } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
153: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
154: } else if (ksp->normtype == KSP_NATURAL_NORM) {
155: dp = sqrt(PetscAbsScalar(beta));
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: i = 0;
165: do {
166: ksp->its = i+1;
167: VecXDot(Z,R,&beta); /* beta <- r'z */
168: if (beta == 0.0) {
169: ksp->reason = KSP_CONVERGED_ATOL;
170: PetscLogInfo(ksp,"KSPSolve_CG:converged due to beta = 0");
171: break;
172: #if !defined(PETSC_USE_COMPLEX)
173: } else if (beta < 0.0) {
174: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
175: PetscLogInfo(ksp,"KSPSolve_CG:diverging due to indefinite preconditioner");
176: break;
177: #endif
178: }
179: if (!i) {
180: VecCopy(Z,P); /* p <- z */
181: } else {
182: b = beta/betaold;
183: if (eigs) {
184: if (ksp->max_it != stored_max_it) {
185: SETERRQ(1,"Can not change maxit AND calculate eigenvalues");
186: }
187: e[i] = sqrt(PetscAbsScalar(b))/a;
188: }
189: VecAYPX(&b,Z,P); /* p <- z + b* p */
190: }
191: betaold = beta;
192: KSP_MatMult(ksp,Amat,P,Z); /* z <- Kp */
193: VecXDot(P,Z,&dpi); /* dpi <- z'p */
194: a = beta/dpi; /* a = beta/p'z */
195: if (eigs) {
196: d[i] = sqrt(PetscAbsScalar(b))*e[i] + 1.0/a;
197: }
198: VecAXPY(&a,P,X); /* x <- x + ap */
199: ma = -a; VecAXPY(&ma,Z,R); /* r <- r - az */
200: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
201: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- Br */
202: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
203: } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
204: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
205: } else if (ksp->normtype == KSP_NATURAL_NORM) {
206: dp = sqrt(PetscAbsScalar(beta));
207: } else {
208: dp = 0.0;
209: }
210: ksp->rnorm = dp;
211: KSPLogResidualHistory(ksp,dp);
212: KSPMonitor(ksp,i+1,dp);
213: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
214: if (ksp->reason) break;
215: if (ksp->normtype != KSP_PRECONDITIONED_NORM) {
216: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- Br */
217: }
218: i++;
219: } while (i<ksp->max_it);
220: if (i == ksp->max_it) {
221: ksp->reason = KSP_DIVERGED_ITS;
222: }
223: return(0);
224: }
225: /*
226: KSPDestroy_CG - Frees all memory space used by the Krylov method
228: */
231: int KSPDestroy_CG(KSP ksp)
232: {
233: KSP_CG *cg = (KSP_CG*)ksp->data;
234: int ierr;
237: /* free space used for singular value calculations */
238: if (ksp->calc_sings) {
239: PetscFree(cg->e);
240: PetscFree(cg->ee);
241: }
243: KSPDefaultFreeWork(ksp);
244:
245: /* free the context variable */
246: PetscFree(cg);
247: return(0);
248: }
250: /*
251: KSPView_CG - Prints information about the current Krylov method being used
253: Currently this only prints information to a file (or stdout) about the
254: symmetry of the problem. If your Krylov method has special options or
255: flags that information should be printed here.
257: */
260: int KSPView_CG(KSP ksp,PetscViewer viewer)
261: {
262: #if defined(PETSC_USE_COMPLEX)
263: KSP_CG *cg = (KSP_CG *)ksp->data;
264: int ierr;
265: PetscTruth isascii;
268: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
269: if (isascii) {
270: if (cg->type == KSP_CG_HERMITIAN) {
271: PetscViewerASCIIPrintf(viewer," CG: variant for complex, Hermitian system\n");
272: } else if (cg->type == KSP_CG_SYMMETRIC) {
273: PetscViewerASCIIPrintf(viewer," CG: variant for complex, symmetric system\n");
274: } else {
275: PetscViewerASCIIPrintf(viewer," CG: unknown variant\n");
276: }
277: } else {
278: SETERRQ1(1,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
279: }
280: #endif
281: return(0);
282: }
284: /*
285: KSPSetFromOptions_CG - Checks the options database for options related to the
286: conjugate gradient method.
287: */
290: int KSPSetFromOptions_CG(KSP ksp)
291: {
292: #if defined(PETSC_USE_COMPLEX)
293: int ierr;
294: PetscTruth flg;
295: #endif
298: #if defined(PETSC_USE_COMPLEX)
299: PetscOptionsHead("KSP CG options");
300: PetscOptionsLogicalGroupBegin("-ksp_cg_Hermitian","Matrix is Hermitian","KSPCGSetType",&flg);
301: if (flg) { KSPCGSetType(ksp,KSP_CG_HERMITIAN); }
302: PetscOptionsLogicalGroupEnd("-ksp_cg_symmetric","Matrix is complex symmetric, not Hermitian","KSPCGSetType",&flg);
303: if (flg) { KSPCGSetType(ksp,KSP_CG_SYMMETRIC); }
304: PetscOptionsTail();
305: #endif
306: return(0);
307: }
309: /*
310: KSPCGSetType_CG - This is an option that is SPECIFIC to this particular Krylov method.
311: This routine is registered below in KSPCreate_CG() and called from the
312: routine KSPCGSetType() (see the file cgtype.c).
314: This must be wrapped in an EXTERN_C_BEGIN to be dynamically linkable in C++
315: */
316: EXTERN_C_BEGIN
319: int KSPCGSetType_CG(KSP ksp,KSPCGType type)
320: {
321: KSP_CG *cg;
324: cg = (KSP_CG *)ksp->data;
325: cg->type = type;
326: return(0);
327: }
328: EXTERN_C_END
330: /*
331: KSPCreate_CG - Creates the data structure for the Krylov method CG and sets the
332: function pointers for all the routines it needs to call (KSPSolve_CG() etc)
334: It must be wrapped in EXTERN_C_BEGIN to be dynamically linkable in C++
335: */
336: /*MC
337: KSPCG - The preconditioned conjugate gradient (PCG) iterative method
339: Options Database Keys:
340: + -ksp_cg_Hermitian - (for complex matrices only) indicates the matrix is Hermitian
341: - -ksp_cg_symmetric - (for complex matrices only) indicates the matrix is symmetric
343: Level: beginner
345: Notes: The PCG method requires both the matrix and preconditioner to
346: be symmetric positive (semi) definite
348: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
349: KSPCGSetType()
351: M*/
352: EXTERN_C_BEGIN
355: int KSPCreate_CG(KSP ksp)
356: {
357: int ierr;
358: KSP_CG *cg;
361: PetscNew(KSP_CG,&cg);
362: PetscMemzero(cg,sizeof(KSP_CG));
363: PetscLogObjectMemory(ksp,sizeof(KSP_CG));
364: #if !defined(PETSC_USE_COMPLEX)
365: cg->type = KSP_CG_SYMMETRIC;
366: #else
367: cg->type = KSP_CG_HERMITIAN;
368: #endif
369: ksp->data = (void*)cg;
370: ksp->pc_side = PC_LEFT;
372: /*
373: Sets the functions that are associated with this data structure
374: (in C++ this is the same as defining virtual functions)
375: */
376: ksp->ops->setup = KSPSetUp_CG;
377: ksp->ops->solve = KSPSolve_CG;
378: ksp->ops->destroy = KSPDestroy_CG;
379: ksp->ops->view = KSPView_CG;
380: ksp->ops->setfromoptions = KSPSetFromOptions_CG;
381: ksp->ops->buildsolution = KSPDefaultBuildSolution;
382: ksp->ops->buildresidual = KSPDefaultBuildResidual;
384: /*
385: Attach the function KSPCGSetType_CG() to this object. The routine
386: KSPCGSetType() checks for this attached function and calls it if it finds
387: it. (Sort of like a dynamic member function that can be added at run time
388: */
389: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG",
390: KSPCGSetType_CG);
391: return(0);
392: }
393: EXTERN_C_END