Actual source code: cgne.c
1: #define PETSCKSP_DLL
3: /*
4: cgimpl.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/cgimpl.h
9: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
10: EXTERN PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);
13: /*
14: KSPSetUp_CGNE - Sets up the workspace needed by the CGNE method.
16: IDENTICAL TO THE CG ONE EXCEPT for one extra work vector!
17: */
20: PetscErrorCode KSPSetUp_CGNE(KSP ksp)
21: {
22: KSP_CG *cgP = (KSP_CG*)ksp->data;
24: PetscInt maxit = ksp->max_it;
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(PETSC_ERR_SUP,"No right preconditioning for KSPCGNE");
33: } else if (ksp->pc_side == PC_SYMMETRIC) {
34: SETERRQ(PETSC_ERR_SUP,"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: PetscMalloc4(maxit+1,PetscScalar,&cgP->e,maxit+1,PetscScalar,&cgP->d,maxit+1,PetscReal,&cgP->ee,maxit+1,PetscReal,&cgP->dd);
47: PetscLogObjectMemory(ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));
48: }
49: return(0);
50: }
52: /*
53: KSPSolve_CGNE - This routine actually applies the conjugate gradient
54: method
56: Input Parameter:
57: . ksp - the Krylov space object that was set to use conjugate gradient, by, for
58: example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);
61: Virtually identical to the KSPSolve_CG, it should definitely reuse the same code.
63: */
66: PetscErrorCode KSPSolve_CGNE(KSP ksp)
67: {
69: PetscInt i,stored_max_it,eigs;
70: PetscScalar dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0;
71: PetscReal dp = 0.0;
72: Vec X,B,Z,R,P,T;
73: KSP_CG *cg;
74: Mat Amat,Pmat;
75: MatStructure pflag;
76: PetscTruth diagonalscale,transpose_pc;
79: PCDiagonalScale(ksp->pc,&diagonalscale);
80: if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
81: PCApplyTransposeExists(ksp->pc,&transpose_pc);
83: cg = (KSP_CG*)ksp->data;
84: eigs = ksp->calc_sings;
85: stored_max_it = ksp->max_it;
86: X = ksp->vec_sol;
87: B = ksp->vec_rhs;
88: R = ksp->work[0];
89: Z = ksp->work[1];
90: P = ksp->work[2];
91: T = ksp->work[3];
93: #if !defined(PETSC_USE_COMPLEX)
94: #define VecXDot(x,y,a) VecDot(x,y,a)
95: #else
96: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
97: #endif
99: if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
100: PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);
102: ksp->its = 0;
103: MatMultTranspose(Amat,B,T);
104: if (!ksp->guess_zero) {
105: KSP_MatMult(ksp,Amat,X,P);
106: KSP_MatMultTranspose(ksp,Amat,P,R);
107: VecAYPX(R,-1.0,T);
108: } else {
109: VecCopy(T,R); /* r <- b (x is 0) */
110: }
111: KSP_PCApply(ksp,R,T);
112: if (transpose_pc) {
113: KSP_PCApplyTranspose(ksp,T,Z);
114: } else {
115: KSP_PCApply(ksp,T,Z);
116: }
118: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
119: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
120: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
121: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
122: } else if (ksp->normtype == KSP_NORM_NATURAL) {
123: VecXDot(Z,R,&beta);
124: dp = sqrt(PetscAbsScalar(beta));
125: } else dp = 0.0;
126: KSPLogResidualHistory(ksp,dp);
127: KSPMonitor(ksp,0,dp); /* call any registered monitor routines */
128: ksp->rnorm = dp;
129: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
130: if (ksp->reason) return(0);
132: i = 0;
133: do {
134: ksp->its = i+1;
135: VecXDot(Z,R,&beta); /* beta <- r'z */
136: if (beta == 0.0) {
137: ksp->reason = KSP_CONVERGED_ATOL;
138: PetscInfo(ksp,"converged due to beta = 0\n");
139: break;
140: #if !defined(PETSC_USE_COMPLEX)
141: } else if (beta < 0.0) {
142: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
143: PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
144: break;
145: #endif
146: }
147: if (!i) {
148: VecCopy(Z,P); /* p <- z */
149: b = 0.0;
150: } else {
151: b = beta/betaold;
152: if (eigs) {
153: if (ksp->max_it != stored_max_it) {
154: SETERRQ(PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
155: }
156: e[i] = sqrt(PetscAbsScalar(b))/a;
157: }
158: VecAYPX(P,b,Z); /* p <- z + b* p */
159: }
160: betaold = beta;
161: MatMult(Amat,P,T);
162: MatMultTranspose(Amat,T,Z);
163: VecXDot(P,Z,&dpi); /* dpi <- z'p */
164: a = beta/dpi; /* a = beta/p'z */
165: if (eigs) {
166: d[i] = sqrt(PetscAbsScalar(b))*e[i] + 1.0/a;
167: }
168: VecAXPY(X,a,P); /* x <- x + ap */
169: VecAXPY(R,-a,Z); /* r <- r - az */
170: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
171: KSP_PCApply(ksp,R,T);
172: if (transpose_pc) {
173: KSP_PCApplyTranspose(ksp,T,Z);
174: } else {
175: KSP_PCApply(ksp,T,Z);
176: }
177: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
178: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
179: VecNorm(R,NORM_2,&dp);
180: } else if (ksp->normtype == KSP_NORM_NATURAL) {
181: dp = sqrt(PetscAbsScalar(beta));
182: } else {
183: dp = 0.0;
184: }
185: ksp->rnorm = dp;
186: KSPLogResidualHistory(ksp,dp);
187: KSPMonitor(ksp,i+1,dp);
188: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
189: if (ksp->reason) break;
190: if (ksp->normtype != KSP_NORM_PRECONDITIONED) {
191: if (transpose_pc) {
192: KSP_PCApplyTranspose(ksp,T,Z);
193: } else {
194: KSP_PCApply(ksp,T,Z);
195: }
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: }
205: /*
206: KSPCreate_CGNE - Creates the data structure for the Krylov method CGNE and sets the
207: function pointers for all the routines it needs to call (KSPSolve_CGNE() etc)
210: */
212: /*MC
213: KSPCGNE - Applies the preconditioned conjugate gradient method to the normal equations
214: without explicitly forming A^t*A
216: Options Database Keys:
217: . -ksp_cg_type <Hermitian or symmetric - (for complex matrices only) indicates the matrix is Hermitian or symmetric
220: Level: beginner
222: Notes: eigenvalue computation routines will return information about the
223: spectrum of A^t*A, rather than A.
225: This is NOT a different algorithm then used with KSPCG, it merely uses that algorithm with the
226: matrix defined by A^t*A and preconditioner defined by B^t*B where B is the preconditioner for A.
228: This method requires that one be apply to apply the transpose of the preconditioner and operator
229: as well as the operator and preconditioner. If the transpose of the preconditioner is not available then
230: the preconditioner is used in its place so one ends up preconditioning A'A with B B. Seems odd?
232: Developer Notes: How is this related to the preconditioned LSQR implementation?
234: This object is subclassed off of KSPCG
236: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
237: KSPCGSetType(), KSPBICG
239: M*/
251: PetscErrorCode KSPCreate_CGNE(KSP ksp)
252: {
254: KSP_CG *cg;
257: PetscNewLog(ksp,KSP_CG,&cg);
258: #if !defined(PETSC_USE_COMPLEX)
259: cg->type = KSP_CG_SYMMETRIC;
260: #else
261: cg->type = KSP_CG_HERMITIAN;
262: #endif
263: ksp->data = (void*)cg;
264: ksp->pc_side = PC_LEFT;
266: /*
267: Sets the functions that are associated with this data structure
268: (in C++ this is the same as defining virtual functions)
269: */
270: ksp->ops->setup = KSPSetUp_CGNE;
271: ksp->ops->solve = KSPSolve_CGNE;
272: ksp->ops->destroy = KSPDestroy_CG;
273: ksp->ops->view = KSPView_CG;
274: ksp->ops->setfromoptions = KSPSetFromOptions_CG;
275: ksp->ops->buildsolution = KSPDefaultBuildSolution;
276: ksp->ops->buildresidual = KSPDefaultBuildResidual;
278: /*
279: Attach the function KSPCGSetType_CGNE() to this object. The routine
280: KSPCGSetType() checks for this attached function and calls it if it finds
281: it. (Sort of like a dynamic member function that can be added at run time
282: */
283: PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPCGSetType_C","KSPCGSetType_CG",KSPCGSetType_CG);
284: return(0);
285: }