2: /*
3: cgimpl.h defines the simple data structured used to store information
4: related to the type of matrix (e.g. complex symmetric) being solved and
5: data used during the optional Lanczo process used to compute eigenvalues
6: */
7: #include <../src/ksp/ksp/impls/cg/cgimpl.h> /*I "petscksp.h" I*/
8: extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal*,PetscReal*);
9: extern PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*);
13: static PetscErrorCode KSPCGSetType_CGNE(KSP ksp,KSPCGType type) 14: {
15: KSP_CG *cg = (KSP_CG*)ksp->data;
18: cg->type = type;
19: return(0);
20: }
23: /*
24: KSPSetUp_CGNE - Sets up the workspace needed by the CGNE method.
26: IDENTICAL TO THE CG ONE EXCEPT for one extra work vector!
27: */
30: PetscErrorCode KSPSetUp_CGNE(KSP ksp) 31: {
32: KSP_CG *cgP = (KSP_CG*)ksp->data;
34: PetscInt maxit = ksp->max_it;
37: /* get work vectors needed by CGNE */
38: KSPSetWorkVecs(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,&cgP->e,maxit+1,&cgP->d,maxit+1,&cgP->ee,maxit+1,&cgP->dd);
47: PetscLogObjectMemory((PetscObject)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: PetscBool diagonalscale,transpose_pc;
78: PCGetDiagonalScale(ksp->pc,&diagonalscale);
79: if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
80: PCApplyTransposeExists(ksp->pc,&transpose_pc);
82: cg = (KSP_CG*)ksp->data;
83: eigs = ksp->calc_sings;
84: stored_max_it = ksp->max_it;
85: X = ksp->vec_sol;
86: B = ksp->vec_rhs;
87: R = ksp->work[0];
88: Z = ksp->work[1];
89: P = ksp->work[2];
90: T = ksp->work[3];
92: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a)) 94: if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
95: PCGetOperators(ksp->pc,&Amat,&Pmat);
97: ksp->its = 0;
98: KSP_MatMultTranspose(ksp,Amat,B,T);
99: if (!ksp->guess_zero) {
100: KSP_MatMult(ksp,Amat,X,P);
101: KSP_MatMultTranspose(ksp,Amat,P,R);
102: VecAYPX(R,-1.0,T);
103: } else {
104: VecCopy(T,R); /* r <- b (x is 0) */
105: }
106: KSP_PCApply(ksp,R,T);
107: if (transpose_pc) {
108: KSP_PCApplyTranspose(ksp,T,Z);
109: } else {
110: KSP_PCApply(ksp,T,Z);
111: }
113: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
114: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
115: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
116: VecNorm(R,NORM_2,&dp); /* dp <- r'*r */
117: } else if (ksp->normtype == KSP_NORM_NATURAL) {
118: VecXDot(Z,R,&beta);
119: dp = PetscSqrtReal(PetscAbsScalar(beta));
120: } else dp = 0.0;
121: KSPLogResidualHistory(ksp,dp);
122: KSPMonitor(ksp,0,dp);
123: ksp->rnorm = dp;
124: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
125: if (ksp->reason) return(0);
127: i = 0;
128: do {
129: ksp->its = i+1;
130: VecXDot(Z,R,&beta); /* beta <- r'z */
131: if (beta == 0.0) {
132: ksp->reason = KSP_CONVERGED_ATOL;
133: PetscInfo(ksp,"converged due to beta = 0\n");
134: break;
135: #if !defined(PETSC_USE_COMPLEX)
136: } else if (beta < 0.0) {
137: ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
138: PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
139: break;
140: #endif
141: }
142: if (!i) {
143: VecCopy(Z,P); /* p <- z */
144: b = 0.0;
145: } else {
146: b = beta/betaold;
147: if (eigs) {
148: if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
149: e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
150: }
151: VecAYPX(P,b,Z); /* p <- z + b* p */
152: }
153: betaold = beta;
154: KSP_MatMult(ksp,Amat,P,T);
155: KSP_MatMultTranspose(ksp,Amat,T,Z);
156: VecXDot(P,Z,&dpi); /* dpi <- z'p */
157: a = beta/dpi; /* a = beta/p'z */
158: if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
159: VecAXPY(X,a,P); /* x <- x + ap */
160: VecAXPY(R,-a,Z); /* r <- r - az */
161: if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
162: KSP_PCApply(ksp,R,T);
163: if (transpose_pc) {
164: KSP_PCApplyTranspose(ksp,T,Z);
165: } else {
166: KSP_PCApply(ksp,T,Z);
167: }
168: VecNorm(Z,NORM_2,&dp); /* dp <- z'*z */
169: } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
170: VecNorm(R,NORM_2,&dp);
171: } else if (ksp->normtype == KSP_NORM_NATURAL) {
172: dp = PetscSqrtReal(PetscAbsScalar(beta));
173: } else {
174: dp = 0.0;
175: }
176: ksp->rnorm = dp;
177: KSPLogResidualHistory(ksp,dp);
178: if (eigs) cg->ned = ksp->its;
179: KSPMonitor(ksp,i+1,dp);
180: (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
181: if (ksp->reason) break;
182: if (ksp->normtype != KSP_NORM_PRECONDITIONED) {
183: if (transpose_pc) {
184: KSP_PCApplyTranspose(ksp,T,Z);
185: } else {
186: KSP_PCApply(ksp,T,Z);
187: }
188: }
189: i++;
190: } while (i<ksp->max_it);
191: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
192: return(0);
193: }
195: /*
196: KSPCreate_CGNE - Creates the data structure for the Krylov method CGNE and sets the
197: function pointers for all the routines it needs to call (KSPSolve_CGNE() etc)
199: It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
200: */
202: /*MC
203: KSPCGNE - Applies the preconditioned conjugate gradient method to the normal equations
204: without explicitly forming A^t*A
206: Options Database Keys:
207: . -ksp_cg_type <Hermitian or symmetric - (for complex matrices only) indicates the matrix is Hermitian or symmetric
210: Level: beginner
212: Notes: eigenvalue computation routines will return information about the
213: spectrum of A^t*A, rather than A.
216: CGNE is a general-purpose non-symmetric method. It works well when the singular values are much better behaved than
217: eigenvalues. A unitary matrix is a classic example where CGNE converges in one iteration, but GMRES and CGS need N
218: iterations (see Nachtigal, Reddy, and Trefethen, "How fast are nonsymmetric matrix iterations", 1992). If you intend
219: to solve least squares problems, use KSPLSQR.
221: This is NOT a different algorithm then used with KSPCG, it merely uses that algorithm with the
222: matrix defined by A^t*A and preconditioner defined by B^t*B where B is the preconditioner for A.
224: This method requires that one be able to apply the transpose of the preconditioner and operator
225: as well as the operator and preconditioner. If the transpose of the preconditioner is not available then
226: the preconditioner is used in its place so one ends up preconditioning A'A with B B. Seems odd?
228: This only supports left preconditioning.
230: This object is subclassed off of KSPCG232: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
233: KSPCGSetType(), KSPBICG235: M*/
237: extern PetscErrorCode KSPDestroy_CG(KSP);
238: extern PetscErrorCode KSPReset_CG(KSP);
239: extern PetscErrorCode KSPView_CG(KSP,PetscViewer);
240: extern PetscErrorCode KSPSetFromOptions_CG(PetscOptionItems *PetscOptionsObject,KSP);
241: PETSC_EXTERN PetscErrorCode KSPCGSetType_CG(KSP,KSPCGType);
245: PETSC_EXTERN PetscErrorCode KSPCreate_CGNE(KSP ksp)246: {
248: KSP_CG *cg;
251: PetscNewLog(ksp,&cg);
252: #if !defined(PETSC_USE_COMPLEX)
253: cg->type = KSP_CG_SYMMETRIC;
254: #else
255: cg->type = KSP_CG_HERMITIAN;
256: #endif
257: ksp->data = (void*)cg;
258: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
259: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
260: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
262: /*
263: Sets the functions that are associated with this data structure
264: (in C++ this is the same as defining virtual functions)
265: */
266: ksp->ops->setup = KSPSetUp_CGNE;
267: ksp->ops->solve = KSPSolve_CGNE;
268: ksp->ops->destroy = KSPDestroy_CG;
269: ksp->ops->view = KSPView_CG;
270: ksp->ops->setfromoptions = KSPSetFromOptions_CG;
271: ksp->ops->buildsolution = KSPBuildSolutionDefault;
272: ksp->ops->buildresidual = KSPBuildResidualDefault;
274: /*
275: Attach the function KSPCGSetType_CGNE() to this object. The routine
276: KSPCGSetType() checks for this attached function and calls it if it finds
277: it. (Sort of like a dynamic member function that can be added at run time
278: */
279: PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",KSPCGSetType_CGNE);
280: return(0);
281: }