Actual source code: cgne.c
1: /*$Id: cgne.c,v 1.117 2001/08/07 03:03:50 balay Exp $*/
3: /*
4: cgctx.h defines the simple data structured used to store information
5: related to the type of matrix (e.g. complex symmetric) being solved and
6: data used during the optional Lanczo process used to compute eigenvalues
7: */
8: #include src/ksp/ksp/impls/cg/cgctx.h
9: EXTERN int KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
10: EXTERN int KSPComputeEigenvalues_CG(KSP,int,PetscReal *,PetscReal *,int *);
13: /*
14: KSPSetUp_CGNE - Sets up the workspace needed by the CGNE method.
16: This is called once, usually automatically by KSPSolve() or KSPSetUp()
17: but can be called directly by KSPSetUp()
18: */
21: int KSPSetUp_CGNE(KSP ksp)
22: {
23: KSP_CG *cgP = (KSP_CG*)ksp->data;
24: int maxit = ksp->max_it,ierr;
27: /*
28: This implementation of CGNE only handles left preconditioning
29: so generate an error otherwise.
30: */
31: if (ksp->pc_side == PC_RIGHT) {
32: SETERRQ(2,"No right preconditioning for KSPCGNE");
33: } else if (ksp->pc_side == PC_SYMMETRIC) {
34: SETERRQ(2,"No symmetric preconditioning for KSPCGNE");
35: }
37: /* get work vectors needed by CGNE */
38: KSPDefaultGetWork(ksp,4);
40: /*
41: If user requested computations of eigenvalues then allocate work
42: work space needed
43: */
44: if (ksp->calc_sings) {
45: /* get space to store tridiagonal matrix for Lanczos */
46: PetscMalloc(2*(maxit+1)*sizeof(PetscScalar),&cgP->e);
47: PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
48: cgP->d = cgP->e + maxit + 1;
49: PetscMalloc(2*(maxit+1)*sizeof(PetscReal),&cgP->ee);
50: PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
51: cgP->dd = cgP->ee + maxit + 1;
52: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
53: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_CG;
54: }
55: return(0);
56: }
58: /*
59: KSPSolve_CGNE - This routine actually applies the conjugate gradient
60: method
62: Input Parameter:
63: . ksp - the Krylov space object that was set to use conjugate gradient, by, for
64: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
66: Output Parameter:
67: . its - number of iterations used
69: */
72: int KSPSolve_CGNE(KSP ksp)
73: {
74: int ierr,i,stored_max_it,eigs;
75: PetscScalar dpi,a = 1.0,beta,betaold = 1.0,b,*e = 0,*d = 0,mone = -1.0,ma;
76: PetscReal dp = 0.0;
77: Vec X,B,Z,R,P,T;
78: KSP_CG *cg;
79: Mat Amat,Pmat;
80: MatStructure pflag;
81: PetscTruth diagonalscale,transpose_pc;
84: PCDiagonalScale(ksp->B,&diagonalscale);
85: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
86: PCHasApplyTranspose(ksp->B,&transpose_pc);
88: cg = (KSP_CG*)ksp->data;
89: eigs = ksp->calc_sings;
90: stored_max_it = ksp->max_it;
91: X = ksp->vec_sol;
92: B = ksp->vec_rhs;
93: R = ksp->work[0];
94: Z = ksp->work[1];
95: P = ksp->work[2];
96: T = ksp->work[3];
98: #if !defined(PETSC_USE_COMPLEX)
99: #define VecXDot(x,y,a) VecDot(x,y,a)
100: #else
101: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
102: #endif
104: if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; b = 0.0; }
105: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
107: ksp->its = 0;
108: MatMultTranspose(Amat,B,T);
109: if (!ksp->guess_zero) {
110: KSP_MatMult(ksp,Amat,X,P);
111: KSP_MatMultTranspose(ksp,Amat,P,R);
112: VecAYPX(&mone,T,R);
113: } else {
114: VecCopy(T,R); /* r <- b (x is 0) */
115: }
116: PCApply(ksp->B,R,T,ksp->pc_side);
117: if (transpose_pc) {
118: PCApplyTranspose(ksp->B,T,Z);
119: } else {
120: PCApply(ksp->B,T,Z,ksp->pc_side);
121: }
123: VecXDot(Z,R,&beta);
124: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
125: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
126: } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
127: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
128: } else if (ksp->normtype == KSP_NATURAL_NORM) {
129: dp = sqrt(PetscAbsScalar(beta));
130: } else dp = 0.0;
131: KSPLogResidualHistory(ksp,dp);
132: KSPMonitor(ksp,0,dp); /* call any registered monitor routines */
133: ksp->rnorm = dp;
134: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
135: if (ksp->reason) return(0);
137: i = 0;
138: do {
139: ksp->its = i+1;
140: VecXDot(Z,R,&beta); /* beta <- r'z */
141: if (beta == 0.0) {
142: ksp->reason = KSP_CONVERGED_ATOL;
143: PetscLogInfo(ksp,"KSPSolve_CGNE:converged due to beta = 0");
144: break;
145: #if !defined(PETSC_USE_COMPLEX)
146: } else if (beta < 0.0) {
147: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
148: PetscLogInfo(ksp,"KSPSolve_CGNE:diverging due to indefinite preconditioner");
149: break;
150: #endif
151: }
152: if (!i) {
153: VecCopy(Z,P); /* p <- z */
154: } else {
155: b = beta/betaold;
156: if (eigs) {
157: if (ksp->max_it != stored_max_it) {
158: SETERRQ(1,"Can not change maxit AND calculate eigenvalues");
159: }
160: e[i] = sqrt(PetscAbsScalar(b))/a;
161: }
162: VecAYPX(&b,Z,P); /* p <- z + b* p */
163: }
164: betaold = beta;
165: MatMult(Amat,P,T);
166: MatMultTranspose(Amat,T,Z);
167: VecXDot(P,Z,&dpi); /* dpi <- z'p */
168: a = beta/dpi; /* a = beta/p'z */
169: if (eigs) {
170: d[i] = sqrt(PetscAbsScalar(b))*e[i] + 1.0/a;
171: }
172: VecAXPY(&a,P,X); /* x <- x + ap */
173: ma = -a; VecAXPY(&ma,Z,R); /* r <- r - az */
174: if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
175: PCApply(ksp->B,R,T,ksp->pc_side);
176: if (transpose_pc) {
177: PCApplyTranspose(ksp->B,T,Z);
178: } else {
179: PCApply(ksp->B,T,Z,ksp->pc_side);
180: }
181: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
182: } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
183: VecNorm(R,NORM_2,&dp);
184: } else if (ksp->normtype == KSP_NATURAL_NORM) {
185: dp = sqrt(PetscAbsScalar(beta));
186: } else {
187: dp = 0.0;
188: }
189: ksp->rnorm = dp;
190: KSPLogResidualHistory(ksp,dp);
191: KSPMonitor(ksp,i+1,dp);
192: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
193: if (ksp->reason) break;
194: if (ksp->normtype != KSP_PRECONDITIONED_NORM) {
195: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- Br */
196: }
197: i++;
198: } while (i<ksp->max_it);
199: if (i == ksp->max_it) {
200: ksp->reason = KSP_DIVERGED_ITS;
201: }
202: return(0);
203: }
204: /*
205: KSPDestroy_CGNE - Frees all memory space used by the Krylov method
207: */
210: int KSPDestroy_CGNE(KSP ksp)
211: {
212: KSP_CG *cg = (KSP_CG*)ksp->data;
213: int ierr;
216: /* free space used for singular value calculations */
217: if (ksp->calc_sings) {
218: PetscFree(cg->e);
219: PetscFree(cg->ee);
220: }
222: KSPDefaultFreeWork(ksp);
223:
224: /* free the context variable */
225: PetscFree(cg);
226: return(0);
227: }
229: /*
230: KSPView_CGNE - Prints information about the current Krylov method being used
232: Currently this only prints information to a file (or stdout) about the
233: symmetry of the problem. If your Krylov method has special options or
234: flags that information should be printed here.
236: */
239: int KSPView_CGNE(KSP ksp,PetscViewer viewer)
240: {
241: #if defined(PETSC_USE_COMPLEX)
242: KSP_CG *cg = (KSP_CG *)ksp->data;
243: int ierr;
244: PetscTruth isascii;
247: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
248: if (isascii) {
249: if (cg->type == KSP_CG_HERMITIAN) {
250: PetscViewerASCIIPrintf(viewer," CG: variant for complex, Hermitian system\n");
251: } else if (cg->type == KSP_CG_SYMMETRIC) {
252: PetscViewerASCIIPrintf(viewer," CG: variant for complex, symmetric system\n");
253: } else {
254: PetscViewerASCIIPrintf(viewer," CG: unknown variant\n");
255: }
256: } else {
257: SETERRQ1(1,"Viewer type %s not supported for KSP cg",((PetscObject)viewer)->type_name);
258: }
259: #endif
260: return(0);
261: }
263: /*
264: KSPSetFromOptions_CGNE - Checks the options database for options related to the
265: conjugate gradient method.
266: */
269: int KSPSetFromOptions_CGNE(KSP ksp)
270: {
271: #if defined(PETSC_USE_COMPLEX)
272: int ierr;
273: PetscTruth flg;
274: #endif
277: #if defined(PETSC_USE_COMPLEX)
278: PetscOptionsHead("KSP CGNE options");
279: PetscOptionsLogicalGroupBegin("-ksp_cgne_Hermitian","Matrix is Hermitian","KSPCGSetType",&flg);
280: if (flg) { KSPCGSetType(ksp,KSP_CG_HERMITIAN); }
281: PetscOptionsLogicalGroupEnd("-ksp_cgne_symmetric","Matrix is complex symmetric, not Hermitian","KSPCGSetType",&flg);
282: if (flg) { KSPCGSetType(ksp,KSP_CG_SYMMETRIC); }
283: PetscOptionsTail();
284: #endif
285: return(0);
286: }
288: /*
289: KSPCGSetType_CGNE - This is an option that is SPECIFIC to this particular Krylov method.
290: This routine is registered below in KSPCreate_CGNE() and called from the
291: routine KSPCGSetType() (see the file cgtype.c).
293: This must be wrapped in an EXTERN_C_BEGIN to be dynamically linkable in C++
294: */
295: EXTERN_C_BEGIN
298: int KSPCGSetType_CGNE(KSP ksp,KSPCGType type)
299: {
300: KSP_CG *cg;
303: cg = (KSP_CG *)ksp->data;
304: cg->type = type;
305: return(0);
306: }
307: EXTERN_C_END
309: /*
310: KSPCreate_CGNE - Creates the data structure for the Krylov method CGNE and sets the
311: function pointers for all the routines it needs to call (KSPSolve_CGNE() etc)
313: It must be wrapped in EXTERN_C_BEGIN to be dynamically linkable in C++
314: */
316: /*MC
317: KSPCGNE - Applies the preconditioned conjugate gradient method to the normal equations
318: without explicitly forming A^t*A
320: Options Database Keys:
321: + -ksp_cgne_Hermitian - (for complex matrices only) indicates the matrix is Hermitian
322: - -ksp_cgne_symmetric - (for complex matrices only) indicates the matrix is symmetric
324: Level: beginner
326: Notes: eigenvalue computation routines will return information about the
327: spectrum of A^tA, rather than A.
329: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
330: KSPCGSetType()
332: M*/
334: EXTERN_C_BEGIN
337: int KSPCreate_CGNE(KSP ksp)
338: {
339: int ierr;
340: KSP_CG *cg;
343: PetscNew(KSP_CG,&cg);
344: PetscMemzero(cg,sizeof(KSP_CG));
345: PetscLogObjectMemory(ksp,sizeof(KSP_CG));
346: #if !defined(PETSC_USE_COMPLEX)
347: cg->type = KSP_CG_SYMMETRIC;
348: #else
349: cg->type = KSP_CG_HERMITIAN;
350: #endif
351: ksp->data = (void*)cg;
352: ksp->pc_side = PC_LEFT;
354: /*
355: Sets the functions that are associated with this data structure
356: (in C++ this is the same as defining virtual functions)
357: */
358: ksp->ops->setup = KSPSetUp_CGNE;
359: ksp->ops->solve = KSPSolve_CGNE;
360: ksp->ops->destroy = KSPDestroy_CGNE;
361: ksp->ops->view = KSPView_CGNE;
362: ksp->ops->setfromoptions = KSPSetFromOptions_CGNE;
363: ksp->ops->buildsolution = KSPDefaultBuildSolution;
364: ksp->ops->buildresidual = KSPDefaultBuildResidual;
366: /*
367: Attach the function KSPCGSetType_CGNE() to this object. The routine
368: KSPCGSetType() checks for this attached function and calls it if it finds
369: it. (Sort of like a dynamic member function that can be added at run time
370: */
371: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CGNE",
372: KSPCGSetType_CGNE);
373: return(0);
374: }
375: EXTERN_C_END