Actual source code: cgne.c

petsc-3.3-p7 2013-05-11
  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 *);


 12: /*
 13:      KSPSetUp_CGNE - Sets up the workspace needed by the CGNE method. 

 15:      IDENTICAL TO THE CG ONE EXCEPT for one extra work vector!
 16: */
 19: PetscErrorCode KSPSetUp_CGNE(KSP ksp)
 20: {
 21:   KSP_CG         *cgP = (KSP_CG*)ksp->data;
 23:   PetscInt       maxit = ksp->max_it;

 26:   /* get work vectors needed by CGNE */
 27:   KSPDefaultGetWork(ksp,4);

 29:   /*
 30:      If user requested computations of eigenvalues then allocate work
 31:      work space needed
 32:   */
 33:   if (ksp->calc_sings) {
 34:     /* get space to store tridiagonal matrix for Lanczos */
 35:     PetscMalloc4(maxit+1,PetscScalar,&cgP->e,maxit+1,PetscScalar,&cgP->d,maxit+1,PetscReal,&cgP->ee,maxit+1,PetscReal,&cgP->dd);
 36:     PetscLogObjectMemory(ksp,2*(maxit+1)*(sizeof(PetscScalar)+sizeof(PetscReal)));
 37:   }
 38:   return(0);
 39: }

 41: /*
 42:        KSPSolve_CGNE - This routine actually applies the conjugate gradient 
 43:     method

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


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

 52: */
 55: PetscErrorCode  KSPSolve_CGNE(KSP ksp)
 56: {
 58:   PetscInt       i,stored_max_it,eigs;
 59:   PetscScalar    dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0;
 60:   PetscReal      dp = 0.0;
 61:   Vec            X,B,Z,R,P,T;
 62:   KSP_CG         *cg;
 63:   Mat            Amat,Pmat;
 64:   MatStructure   pflag;
 65:   PetscBool      diagonalscale,transpose_pc;

 68:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 69:   if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);
 70:   PCApplyTransposeExists(ksp->pc,&transpose_pc);

 72:   cg            = (KSP_CG*)ksp->data;
 73:   eigs          = ksp->calc_sings;
 74:   stored_max_it = ksp->max_it;
 75:   X             = ksp->vec_sol;
 76:   B             = ksp->vec_rhs;
 77:   R             = ksp->work[0];
 78:   Z             = ksp->work[1];
 79:   P             = ksp->work[2];
 80:   T             = ksp->work[3];

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

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

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

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

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

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

194:     It must be wrapped in EXTERN_C_BEGIN to be dynamically linkable in C++
195: */

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

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


205:    Level: beginner

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

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

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

217:    This only supports left preconditioning.

219:    Developer Notes: How is this related to the preconditioned LSQR implementation?

221:    This object is subclassed off of KSPCG

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

226: M*/

228: extern PetscErrorCode KSPDestroy_CG(KSP);
229: extern PetscErrorCode KSPReset_CG(KSP);
230: extern PetscErrorCode KSPView_CG(KSP,PetscViewer);
231: extern PetscErrorCode KSPSetFromOptions_CG(KSP);
232: EXTERN_C_BEGIN
233: extern PetscErrorCode  KSPCGSetType_CG(KSP,KSPCGType);
234: EXTERN_C_END

236: EXTERN_C_BEGIN
239: PetscErrorCode  KSPCreate_CGNE(KSP ksp)
240: {
242:   KSP_CG         *cg;

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

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

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