Actual source code: cgne.c


  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>
  8: extern PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal*,PetscReal*);
  9: extern PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal*,PetscReal*,PetscInt*);

 11: static PetscErrorCode  KSPCGSetType_CGNE(KSP ksp,KSPCGType type)
 12: {
 13:   KSP_CG *cg = (KSP_CG*)ksp->data;

 16:   cg->type = type;
 17:   return(0);
 18: }

 20: /*
 21:      KSPSetUp_CGNE - Sets up the workspace needed by the CGNE method.

 23:      IDENTICAL TO THE CG ONE EXCEPT for one extra work vector!
 24: */
 25: static PetscErrorCode KSPSetUp_CGNE(KSP ksp)
 26: {
 27:   KSP_CG         *cgP = (KSP_CG*)ksp->data;
 29:   PetscInt       maxit = ksp->max_it;

 32:   /* get work vectors needed by CGNE */
 33:   KSPSetWorkVecs(ksp,4);

 35:   /*
 36:      If user requested computations of eigenvalues then allocate work
 37:      work space needed
 38:   */
 39:   if (ksp->calc_sings) {
 40:     /* get space to store tridiagonal matrix for Lanczos */
 41:     PetscMalloc4(maxit+1,&cgP->e,maxit+1,&cgP->d,maxit+1,&cgP->ee,maxit+1,&cgP->dd);
 42:     PetscLogObjectMemory((PetscObject)ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));

 44:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 45:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 46:   }
 47:   return(0);
 48: }

 50: /*
 51:        KSPSolve_CGNE - This routine actually applies the conjugate gradient
 52:     method

 54:    Input Parameter:
 55: .     ksp - the Krylov space object that was set to use conjugate gradient, by, for
 56:             example, KSPCreate(MPI_Comm,KSP *ksp); KSPSetType(ksp,KSPCG);

 58:     Virtually identical to the KSPSolve_CG, it should definitely reuse the same code.

 60: */
 61: static PetscErrorCode  KSPSolve_CGNE(KSP ksp)
 62: {
 64:   PetscInt       i,stored_max_it,eigs;
 65:   PetscScalar    dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = NULL,*d = NULL;
 66:   PetscReal      dp = 0.0;
 67:   Vec            X,B,Z,R,P,T;
 68:   KSP_CG         *cg;
 69:   Mat            Amat,Pmat;
 70:   PetscBool      diagonalscale,transpose_pc;

 73:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 74:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
 75:   PCApplyTransposeExists(ksp->pc,&transpose_pc);

 77:   cg            = (KSP_CG*)ksp->data;
 78:   eigs          = ksp->calc_sings;
 79:   stored_max_it = ksp->max_it;
 80:   X             = ksp->vec_sol;
 81:   B             = ksp->vec_rhs;
 82:   R             = ksp->work[0];
 83:   Z             = ksp->work[1];
 84:   P             = ksp->work[2];
 85:   T             = ksp->work[3];

 87: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))

 89:   if (eigs) {e = cg->e; d = cg->d; e[0] = 0.0; }
 90:   PCGetOperators(ksp->pc,&Amat,&Pmat);

 92:   ksp->its = 0;
 93:   KSP_MatMultTranspose(ksp,Amat,B,T);
 94:   if (!ksp->guess_zero) {
 95:     KSP_MatMult(ksp,Amat,X,P);
 96:     KSP_MatMultTranspose(ksp,Amat,P,R);
 97:     VecAYPX(R,-1.0,T);
 98:   } else {
 99:     VecCopy(T,R);              /*     r <- b (x is 0) */
100:   }
101:   if (transpose_pc) {
102:     KSP_PCApplyTranspose(ksp,R,T);
103:   } else {
104:     KSP_PCApply(ksp,R,T);
105:   }
106:   KSP_PCApply(ksp,T,Z);

108:   if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
109:     VecNorm(Z,NORM_2,&dp); /*    dp <- z'*z       */
110:   } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
111:     VecNorm(R,NORM_2,&dp); /*    dp <- r'*r       */
112:   } else if (ksp->normtype == KSP_NORM_NATURAL) {
113:     VecXDot(Z,R,&beta);
114:     KSPCheckDot(ksp,beta);
115:     dp   = PetscSqrtReal(PetscAbsScalar(beta));
116:   } else dp = 0.0;
117:   KSPLogResidualHistory(ksp,dp);
118:   KSPMonitor(ksp,0,dp);
119:   ksp->rnorm = dp;
120:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
121:   if (ksp->reason) return(0);

123:   i = 0;
124:   do {
125:     ksp->its = i+1;
126:     VecXDot(Z,R,&beta); /*     beta <- r'z     */
127:     KSPCheckDot(ksp,beta);
128:     if (beta == 0.0) {
129:       ksp->reason = KSP_CONVERGED_ATOL;
130:       PetscInfo(ksp,"converged due to beta = 0\n");
131:       break;
132: #if !defined(PETSC_USE_COMPLEX)
133:     } else if (beta < 0.0) {
134:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
135:       PetscInfo(ksp,"diverging due to indefinite preconditioner\n");
136:       break;
137: #endif
138:     }
139:     if (!i) {
140:       VecCopy(Z,P);          /*     p <- z          */
141:       b    = 0.0;
142:     } else {
143:       b = beta/betaold;
144:       if (eigs) {
145:         if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues");
146:         e[i] = PetscSqrtReal(PetscAbsScalar(b))/a;
147:       }
148:       VecAYPX(P,b,Z);     /*     p <- z + b* p   */
149:     }
150:     betaold = beta;
151:     KSP_MatMult(ksp,Amat,P,T);
152:     KSP_MatMultTranspose(ksp,Amat,T,Z);
153:     VecXDot(P,Z,&dpi);    /*     dpi <- z'p      */
154:     KSPCheckDot(ksp,dpi);
155:     a       = beta/dpi;                            /*     a = beta/p'z    */
156:     if (eigs) d[i] = PetscSqrtReal(PetscAbsScalar(b))*e[i] + 1.0/a;
157:     VecAXPY(X,a,P);           /*     x <- x + ap     */
158:     VecAXPY(R,-a,Z);                       /*     r <- r - az     */
159:     if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
160:       if (transpose_pc) {
161:         KSP_PCApplyTranspose(ksp,R,T);
162:       } else {
163:         KSP_PCApply(ksp,R,T);
164:       }
165:       KSP_PCApply(ksp,T,Z);
166:       VecNorm(Z,NORM_2,&dp);              /*    dp <- z'*z       */
167:     } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
168:       VecNorm(R,NORM_2,&dp);
169:     } else if (ksp->normtype == KSP_NORM_NATURAL) {
170:       dp = PetscSqrtReal(PetscAbsScalar(beta));
171:     } else dp = 0.0;
172:     ksp->rnorm = dp;
173:     KSPLogResidualHistory(ksp,dp);
174:     if (eigs) cg->ned = ksp->its;
175:     KSPMonitor(ksp,i+1,dp);
176:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
177:     if (ksp->reason) break;
178:     if (ksp->normtype != KSP_NORM_PRECONDITIONED) {
179:       if (transpose_pc) {
180:         KSP_PCApplyTranspose(ksp,R,T);
181:       } else {
182:         KSP_PCApply(ksp,R,T);
183:       }
184:       KSP_PCApply(ksp,T,Z);
185:     }
186:     i++;
187:   } while (i<ksp->max_it);
188:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
189:   return(0);
190: }

192: /*
193:     KSPCreate_CGNE - Creates the data structure for the Krylov method CGNE and sets the
194:        function pointers for all the routines it needs to call (KSPSolve_CGNE() etc)

196:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
197: */

199: /*MC
200:      KSPCGNE - Applies the preconditioned conjugate gradient method to the normal equations
201:           without explicitly forming A^t*A

203:    Options Database Keys:
204: .   -ksp_cg_type <Hermitian or symmetric - (for complex matrices only) indicates the matrix is Hermitian or symmetric

206:    Level: beginner

208:    Notes:
209:     eigenvalue computation routines will return information about the
210:           spectrum of A^t*A, rather than A.

212:    CGNE is a general-purpose non-symmetric method. It works well when the singular values are much better behaved than
213:    eigenvalues. A unitary matrix is a classic example where CGNE converges in one iteration, but GMRES and CGS need N
214:    iterations (see Nachtigal, Reddy, and Trefethen, "How fast are nonsymmetric matrix iterations", 1992). If you intend
215:    to solve least squares problems, use KSPLSQR.

217:    This is NOT a different algorithm than used with KSPCG, it merely uses that algorithm with the
218:    matrix defined by A^t*A and preconditioner defined by B^t*B where B is the preconditioner for A.

220:    This method requires that one be able to apply the transpose of the preconditioner and operator
221:    as well as the operator and preconditioner. If the transpose of the preconditioner is not available then
222:    the preconditioner is used in its place so one ends up preconditioning A'A with B B. Seems odd?

224:    This only supports left preconditioning.

226:    This object is subclassed off of KSPCG

228: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
229:            KSPCGSetType(), KSPBICG

231: M*/

233: PETSC_EXTERN PetscErrorCode KSPCreate_CGNE(KSP ksp)
234: {
236:   KSP_CG         *cg;

239:   PetscNewLog(ksp,&cg);
240: #if !defined(PETSC_USE_COMPLEX)
241:   cg->type = KSP_CG_SYMMETRIC;
242: #else
243:   cg->type = KSP_CG_HERMITIAN;
244: #endif
245:   ksp->data = (void*)cg;
246:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
247:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
248:   KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
249:   KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);

251:   /*
252:        Sets the functions that are associated with this data structure
253:        (in C++ this is the same as defining virtual functions)
254:   */
255:   ksp->ops->setup          = KSPSetUp_CGNE;
256:   ksp->ops->solve          = KSPSolve_CGNE;
257:   ksp->ops->destroy        = KSPDestroy_CG;
258:   ksp->ops->view           = KSPView_CG;
259:   ksp->ops->setfromoptions = KSPSetFromOptions_CG;
260:   ksp->ops->buildsolution  = KSPBuildSolutionDefault;
261:   ksp->ops->buildresidual  = KSPBuildResidualDefault;

263:   /*
264:       Attach the function KSPCGSetType_CGNE() to this object. The routine
265:       KSPCGSetType() checks for this attached function and calls it if it finds
266:       it. (Sort of like a dynamic member function that can be added at run time
267:   */
268:   PetscObjectComposeFunction((PetscObject)ksp,"KSPCGSetType_C",KSPCGSetType_CGNE);
269:   return(0);
270: }