Actual source code: cgne.c

  1: #define PETSCKSP_DLL

  3: /*
  4:        cgctx.h defines the simple data structured used to store information
  5:     related to the type of matrix (e.g. complex symmetric) being solved and
  6:     data used during the optional Lanczo process used to compute eigenvalues
  7: */
 8:  #include src/ksp/ksp/impls/cg/cgctx.h
  9: EXTERN PetscErrorCode KSPComputeExtremeSingularValues_CG(KSP,PetscReal *,PetscReal *);
 10: EXTERN PetscErrorCode KSPComputeEigenvalues_CG(KSP,PetscInt,PetscReal *,PetscReal *,PetscInt *);


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

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

 27:   /* 
 28:        This implementation of CGNE only handles left preconditioning
 29:      so generate an error otherwise.
 30:   */
 31:   if (ksp->pc_side == PC_RIGHT) {
 32:     SETERRQ(PETSC_ERR_SUP,"No right preconditioning for KSPCGNE");
 33:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 34:     SETERRQ(PETSC_ERR_SUP,"No symmetric preconditioning for KSPCGNE");
 35:   }

 37:   /* get work vectors needed by CGNE */
 38:   KSPDefaultGetWork(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:     PetscMalloc(2*(maxit+1)*sizeof(PetscScalar),&cgP->e);
 47:     PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
 48:     cgP->d                         = cgP->e + maxit + 1;
 49:     PetscMalloc(2*(maxit+1)*sizeof(PetscReal),&cgP->ee);
 50:     PetscLogObjectMemory(ksp,2*(maxit+1)*sizeof(PetscScalar));
 51:     cgP->dd                        = cgP->ee + maxit + 1;
 52:     ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_CG;
 53:     ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_CG;
 54:   }
 55:   return(0);
 56: }

 58: /*
 59:        KSPSolve_CGNE - This routine actually applies the conjugate gradient 
 60:     method

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

 66:     Probably virtually identical to the KSPSolve_CG, would be nice if we could reuse the code

 68: */
 71: PetscErrorCode  KSPSolve_CGNE(KSP ksp)
 72: {
 74:   PetscInt       i,stored_max_it,eigs;
 75:   PetscScalar    dpi,a = 1.0,beta,betaold = 1.0,b = 0,*e = 0,*d = 0,mone = -1.0,ma;
 76:   PetscReal      dp = 0.0;
 77:   Vec            X,B,Z,R,P,T;
 78:   KSP_CG         *cg;
 79:   Mat            Amat,Pmat;
 80:   MatStructure   pflag;
 81:   PetscTruth     diagonalscale,transpose_pc;

 84:   PCDiagonalScale(ksp->pc,&diagonalscale);
 85:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);
 86:   PCHasApplyTranspose(ksp->pc,&transpose_pc);

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

 98: #if !defined(PETSC_USE_COMPLEX)
 99: #define VecXDot(x,y,a) VecDot(x,y,a)
100: #else
101: #define VecXDot(x,y,a) (((cg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a))
102: #endif

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

107:   ksp->its = 0;
108:   MatMultTranspose(Amat,B,T);
109:   if (!ksp->guess_zero) {
110:     KSP_MatMult(ksp,Amat,X,P);
111:     KSP_MatMultTranspose(ksp,Amat,P,R);
112:     VecAYPX(R,mone,T);
113:   } else {
114:     VecCopy(T,R);              /*     r <- b (x is 0) */
115:   }
116:   KSP_PCApply(ksp,R,T);
117:   if (transpose_pc) {
118:     KSP_PCApplyTranspose(ksp,T,Z);
119:   } else {
120:     KSP_PCApply(ksp,T,Z);
121:   }

123:   VecXDot(Z,R,&beta);
124:   if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
125:     VecNorm(Z,NORM_2,&dp); /*    dp <- z'*z       */
126:   } else if (ksp->normtype == KSP_UNPRECONDITIONED_NORM) {
127:     VecNorm(R,NORM_2,&dp); /*    dp <- r'*r       */
128:   } else if (ksp->normtype == KSP_NATURAL_NORM) {
129:     dp = sqrt(PetscAbsScalar(beta));
130:   } else dp = 0.0;
131:   KSPLogResidualHistory(ksp,dp);
132:   KSPMonitor(ksp,0,dp);                              /* call any registered monitor routines */
133:   ksp->rnorm = dp;
134:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);      /* test for convergence */
135:   if (ksp->reason) return(0);

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

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

210: */

212: /*MC
213:      KSPCGNE - Applies the preconditioned conjugate gradient method to the normal equations
214:           without explicitly forming A^t*A

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


220:    Level: beginner

222:    Notes: eigenvalue computation routines will return information about the
223:    spectrum of A^tA, rather than A.

225:           This object is subclassed off of KSPCG

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

230: M*/


242: PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate_CGNE(KSP ksp)
243: {
245:   KSP_CG         *cg;

248:   PetscNew(KSP_CG,&cg);
249:   PetscLogObjectMemory(ksp,sizeof(KSP_CG));
250: #if !defined(PETSC_USE_COMPLEX)
251:   cg->type                       = KSP_CG_SYMMETRIC;
252: #else
253:   cg->type                       = KSP_CG_HERMITIAN;
254: #endif
255:   ksp->data                      = (void*)cg;
256:   ksp->pc_side                   = PC_LEFT;

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

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