Actual source code: cgne.c

petsc-master 2016-07-27
Report Typos and Errors
  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: }

 22: /*
 23:      KSPSetUp_CGNE - Sets up the workspace needed by the CGNE method.

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

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

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

 48:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 49:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 50:   }
 51:   return(0);
 52: }

 54: /*
 55:        KSPSolve_CGNE - This routine actually applies the conjugate gradient
 56:     method

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


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

 65: */
 68: static PetscErrorCode  KSPSolve_CGNE(KSP ksp)
 69: {
 71:   PetscInt       i,stored_max_it,eigs;
 72:   PetscScalar    dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0;
 73:   PetscReal      dp = 0.0;
 74:   Vec            X,B,Z,R,P,T;
 75:   KSP_CG         *cg;
 76:   Mat            Amat,Pmat;
 77:   PetscBool      diagonalscale,transpose_pc;

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

 84:   cg            = (KSP_CG*)ksp->data;
 85:   eigs          = ksp->calc_sings;
 86:   stored_max_it = ksp->max_it;
 87:   X             = ksp->vec_sol;
 88:   B             = ksp->vec_rhs;
 89:   R             = ksp->work[0];
 90:   Z             = ksp->work[1];
 91:   P             = ksp->work[2];
 92:   T             = ksp->work[3];

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

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

 99:   ksp->its = 0;
100:   KSP_MatMultTranspose(ksp,Amat,B,T);
101:   if (!ksp->guess_zero) {
102:     KSP_MatMult(ksp,Amat,X,P);
103:     KSP_MatMultTranspose(ksp,Amat,P,R);
104:     VecAYPX(R,-1.0,T);
105:   } else {
106:     VecCopy(T,R);              /*     r <- b (x is 0) */
107:   }
108:   if (transpose_pc) {
109:     KSP_PCApplyTranspose(ksp,R,T);
110:   } else {
111:     KSP_PCApply(ksp,R,T);
112:   }
113:   KSP_PCApply(ksp,T,Z);

115:   if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
116:     VecNorm(Z,NORM_2,&dp); /*    dp <- z'*z       */
117:   } else if (ksp->normtype == KSP_NORM_UNPRECONDITIONED) {
118:     VecNorm(R,NORM_2,&dp); /*    dp <- r'*r       */
119:   } else if (ksp->normtype == KSP_NORM_NATURAL) {
120:     VecXDot(Z,R,&beta);
121:     KSPCheckDot(ksp,beta);
122:     dp   = PetscSqrtReal(PetscAbsScalar(beta));
123:   } else dp = 0.0;
124:   KSPLogResidualHistory(ksp,dp);
125:   KSPMonitor(ksp,0,dp);
126:   ksp->rnorm = dp;
127:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP); /* test for convergence */
128:   if (ksp->reason) return(0);

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

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

205:     It must be labeled as PETSC_EXTERN to be dynamically linkable in C++
206: */

208: /*MC
209:      KSPCGNE - Applies the preconditioned conjugate gradient method to the normal equations
210:           without explicitly forming A^t*A

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


216:    Level: beginner

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


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

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

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

234:    This only supports left preconditioning.

236:    This object is subclassed off of KSPCG

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

241: M*/

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: }