Actual source code: cgs.c

  1: /*                       

  3:     Note that for the complex numbers version, the VecDot() arguments
  4:     within the code MUST remain in the order given for correct computation
  5:     of inner products.
  6: */
 7:  #include src/ksp/ksp/kspimpl.h

 11: static PetscErrorCode KSPSetUp_CGS(KSP ksp)
 12: {

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

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

 33:   PCDiagonalScale(ksp->pc,&diagonalscale);
 34:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);

 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:   VecNorm(R,NORM_2,&dp);
 52:   if (ksp->normtype == KSP_NATURAL_NORM) {
 53:     dp *= dp;
 54:   }
 55:   PetscObjectTakeAccess(ksp);
 56:   ksp->its   = 0;
 57:   ksp->rnorm = dp;
 58:   PetscObjectGrantAccess(ksp);
 59:   KSPLogResidualHistory(ksp,dp);
 60:   KSPMonitor(ksp,0,dp);
 61:   (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
 62:   if (ksp->reason) return(0);

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

 85:   /* Set the initial conditions */
 86:   VecDot(R,RP,&rhoold);        /* rhoold = (r,rp)      */
 87:   VecCopy(R,U);
 88:   VecCopy(R,P);
 89:   KSP_PCApplyBAorAB(ksp,P,V,T);

 91:   i = 0;
 92:   do {

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

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

118:     b    = rho / rhoold;                             /* b <- rho / rhoold    */
119:     VecWAXPY(&b,Q,R,U);        /* u <- r + b q         */
120:     VecAXPY(&b,P,Q);
121:     VecWAXPY(&b,Q,U,P);        /* p <- u + b(q + b p)  */
122:     KSP_PCApplyBAorAB(ksp,P,V,Q);      /* v <- K p    */
123:     rhoold = rho;
124:     i++;
125:   } while (i<ksp->max_it);
126:   if (i >= ksp->max_it) {
127:     ksp->reason = KSP_DIVERGED_ITS;
128:   }

130:   KSPUnwindPreconditioner(ksp,X,T);
131:   return(0);
132: }

134: /*MC
135:      KSPCGS - This code implements the CGS (Conjugate Gradient Squared) method. 

137:    Options Database Keys:
138: .   see KSPSolve()

140:    Level: beginner

142:    Notes: Reference: Sonneveld, 1989.

144: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBCGS
145: M*/
149: PetscErrorCode KSPCreate_CGS(KSP ksp)
150: {
152:   ksp->data                      = (void*)0;
153:   ksp->pc_side                   = PC_LEFT;
154:   ksp->ops->setup                = KSPSetUp_CGS;
155:   ksp->ops->solve                = KSPSolve_CGS;
156:   ksp->ops->destroy              = KSPDefaultDestroy;
157:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
158:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
159:   ksp->ops->setfromoptions       = 0;
160:   ksp->ops->view                 = 0;
161:   return(0);
162: }