Actual source code: cgne.c

petsc-3.9.4 2018-09-11
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>
  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);


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

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

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

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

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

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

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

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

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

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

235: M*/

237: PETSC_EXTERN PetscErrorCode KSPCreate_CGNE(KSP ksp)
238: {
240:   KSP_CG         *cg;

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

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

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