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,8);
 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,BINVF,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:   BINVF   = ksp->work[6];
 45:   U       = ksp->work[7];
 46:   AUQ     = V;

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

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

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

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

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

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

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

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

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

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