Actual source code: cgs.c

  1: /*$Id: cgs.c,v 1.64 2001/08/07 03:03:51 balay Exp $*/

  3: /*                       
  4:     This code implements the CGS (Conjugate Gradient Squared) method. 
  5:     Reference: Sonneveld, 1989.

  7:     Note that for the complex numbers version, the VecDot() arguments
  8:     within the code MUST remain in the order given for correct computation
  9:     of inner products.
 10: */
 11:  #include src/sles/ksp/kspimpl.h

 13: static int KSPSetUp_CGS(KSP ksp)
 14: {

 18:   if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(2,"no symmetric preconditioning for KSPCGS");
 19:   KSPDefaultGetWork(ksp,7);
 20:   return(0);
 21: }

 23: static int  KSPSolve_CGS(KSP ksp,int *its)
 24: {
 25:   int          i,maxit,ierr;
 26:   PetscScalar  rho,rhoold,a,s,b,tmp,one = 1.0;
 27:   Vec          X,B,V,P,R,RP,T,Q,U,AUQ;
 28:   PetscReal    dp = 0.0;
 29:   PetscTruth   diagonalscale;

 32:   ierr    = PCDiagonalScale(ksp->B,&diagonalscale);
 33:   if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);

 35:   maxit   = ksp->max_it;
 36:   X       = ksp->vec_sol;
 37:   B       = ksp->vec_rhs;
 38:   R       = ksp->work[0];
 39:   RP      = ksp->work[1];
 40:   V       = ksp->work[2];
 41:   T       = ksp->work[3];
 42:   Q       = ksp->work[4];
 43:   P       = ksp->work[5];
 44:   U       = ksp->work[6];
 45:   AUQ     = V;

 47:   /* Compute initial preconditioned residual */
 48:   KSPInitialResidual(ksp,X,V,T,R,B);

 50:   /* Test for nothing to do */
 51:   if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
 52:     VecNorm(R,NORM_2,&dp);
 53:   } else if (ksp->normtype == KSP_NATURAL_NORM) {
 54:     VecNorm(R,NORM_2,&dp);
 55:   }
 56:   PetscObjectTakeAccess(ksp);
 57:   ksp->its   = 0;
 58:   ksp->rnorm = dp;
 59:   PetscObjectGrantAccess(ksp);
 60:   KSPLogResidualHistory(ksp,dp);
 61:   KSPMonitor(ksp,0,dp);
 62:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 63:   if (ksp->reason) {*its = 0; return(0);}

 65:   /* Make the initial Rp == R */
 66:   VecCopy(R,RP);

 68:   /* Set the initial conditions */
 69:   VecDot(R,RP,&rhoold);        /* rhoold = (r,rp)      */
 70:   VecCopy(R,U);
 71:   VecCopy(R,P);
 72:   KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T);

 74:   for (i=0; i<maxit; i++) {

 76:     VecDot(V,RP,&s);           /* s <- (v,rp)          */
 77:     a = rhoold / s;                                  /* a <- rho / s         */
 78:     tmp = -a;
 79:     VecWAXPY(&tmp,V,U,Q);      /* q <- u - a v         */
 80:     VecWAXPY(&one,U,Q,T);      /* t <- u + q           */
 81:     VecAXPY(&a,T,X);           /* x <- x + a (u + q)   */
 82:     KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,T,AUQ,U);
 83:     VecAXPY(&tmp,AUQ,R);       /* r <- r - a K (u + q) */
 84:     VecDot(R,RP,&rho);         /* rho <- (r,rp)        */
 85:     if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
 86:       VecNorm(R,NORM_2,&dp);
 87:     } else if (ksp->normtype == KSP_NATURAL_NORM) {
 88:       dp = PetscAbsScalar(rho);
 89:     }

 91:     PetscObjectTakeAccess(ksp);
 92:     ksp->its++;
 93:     ksp->rnorm = dp;
 94:     PetscObjectGrantAccess(ksp);
 95:     KSPLogResidualHistory(ksp,dp);
 96:     KSPMonitor(ksp,i+1,dp);
 97:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
 98:     if (ksp->reason) break;

100:     b    = rho / rhoold;                             /* b <- rho / rhoold    */
101:     VecWAXPY(&b,Q,R,U);        /* u <- r + b q         */
102:     VecAXPY(&b,P,Q);
103:     VecWAXPY(&b,Q,U,P);        /* p <- u + b(q + b p)  */
104:     KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,Q);      /* v <- K p    */
105:     rhoold = rho;
106:   }
107:   if (i == maxit) {
108:     i--;
109:     ksp->reason = KSP_DIVERGED_ITS;
110:   }
111:   *its = i+1;

113:   KSPUnwindPreconditioner(ksp,X,T);
114:   return(0);
115: }

117: EXTERN_C_BEGIN
118: int KSPCreate_CGS(KSP ksp)
119: {
121:   ksp->data                      = (void*)0;
122:   ksp->pc_side                   = PC_LEFT;
123:   ksp->ops->setup                = KSPSetUp_CGS;
124:   ksp->ops->solve                = KSPSolve_CGS;
125:   ksp->ops->destroy              = KSPDefaultDestroy;
126:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
127:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
128:   ksp->ops->setfromoptions       = 0;
129:   ksp->ops->view                 = 0;
130:   return(0);
131: }
132: EXTERN_C_END