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: #undef __FUNCT__  
 15: static int KSPSetUp_CGS(KSP ksp)
 16: {

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

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

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

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

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

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

 68:   /* Make the initial Rp == R */
 69:   VecCopy(R,RP);
 70:   /*  added for Fidap */
 71:   /* Penalize Startup - Isaac Hasbani Trick for CGS 
 72:      Since most initial conditions result in a mostly 0 residual,
 73:      we change all the 0 values in the vector RP to the maximum.
 74:   */
 75:   if (ksp->normtype == KSP_NATURAL_NORM) {
 76:      PetscReal vr0max;
 77:      PetscScalar *tmp_RP=0;
 78:      int         numnp=0, *max_pos=0;
 79:      VecMax(RP, max_pos, &vr0max);
 80:      VecGetArray(RP, &tmp_RP);
 81:      VecGetLocalSize(RP, &numnp);
 82:      for (i=0; i<numnp; i++) {
 83:        if (tmp_RP[i] == 0.0) tmp_RP[i] = vr0max;
 84:      }
 85:      VecRestoreArray(RP, &tmp_RP);
 86:   }
 87:   /*  end of addition for Fidap */

 89:   /* Set the initial conditions */
 90:   VecDot(R,RP,&rhoold);        /* rhoold = (r,rp)      */
 91:   VecCopy(R,U);
 92:   VecCopy(R,P);
 93:   KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T);

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

 97:     VecDot(V,RP,&s);           /* s <- (v,rp)          */
 98:     a = rhoold / s;                                  /* a <- rho / s         */
 99:     tmp = -a;
100:     VecWAXPY(&tmp,V,U,Q);      /* q <- u - a v         */
101:     VecWAXPY(&one,U,Q,T);      /* t <- u + q           */
102:     VecAXPY(&a,T,X);           /* x <- x + a (u + q)   */
103:     KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,T,AUQ,U);
104:     VecAXPY(&tmp,AUQ,R);       /* r <- r - a K (u + q) */
105:     VecDot(R,RP,&rho);         /* rho <- (r,rp)        */
106:     if (ksp->normtype == KSP_PRECONDITIONED_NORM) {
107:       VecNorm(R,NORM_2,&dp);
108:     } else if (ksp->normtype == KSP_NATURAL_NORM) {
109:       dp = PetscAbsScalar(rho);
110:     }

112:     PetscObjectTakeAccess(ksp);
113:     ksp->its++;
114:     ksp->rnorm = dp;
115:     PetscObjectGrantAccess(ksp);
116:     KSPLogResidualHistory(ksp,dp);
117:     KSPMonitor(ksp,i+1,dp);
118:     (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);
119:     if (ksp->reason) break;

121:     b    = rho / rhoold;                             /* b <- rho / rhoold    */
122:     VecWAXPY(&b,Q,R,U);        /* u <- r + b q         */
123:     VecAXPY(&b,P,Q);
124:     VecWAXPY(&b,Q,U,P);        /* p <- u + b(q + b p)  */
125:     KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,Q);      /* v <- K p    */
126:     rhoold = rho;
127:   }
128:   if (i == maxit) {
129:     ksp->reason = KSP_DIVERGED_ITS;
130:   }
131:   *its = ksp->its;

133:   KSPUnwindPreconditioner(ksp,X,T);
134:   return(0);
135: }

137: EXTERN_C_BEGIN
138: #undef __FUNCT__  
140: int KSPCreate_CGS(KSP ksp)
141: {
143:   ksp->data                      = (void*)0;
144:   ksp->pc_side                   = PC_LEFT;
145:   ksp->ops->setup                = KSPSetUp_CGS;
146:   ksp->ops->solve                = KSPSolve_CGS;
147:   ksp->ops->destroy              = KSPDefaultDestroy;
148:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
149:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
150:   ksp->ops->setfromoptions       = 0;
151:   ksp->ops->view                 = 0;
152:   return(0);
153: }
154: EXTERN_C_END