Actual source code: symmlq.c

  1: /*$Id: symmlq.c,v 1.16 2001/08/07 03:03:56 balay Exp $*/
  2: /*                       
  3:     This code implements the SYMMLQ method. 
  4:     Reference: Paige & Saunders, 1975.

  6:     Contributed by: Hong Zhang
  7: */
 8:  #include src/sles/ksp/kspimpl.h

 10: typedef struct {
 11:   PetscReal haptol;
 12: } KSP_SYMMLQ;

 14: int KSPSetUp_SYMMLQ(KSP ksp)
 15: {

 19:   if (ksp->pc_side == PC_RIGHT) {
 20:     SETERRQ(2,"No right preconditioning for KSPSYMMLQ");
 21:   } else if (ksp->pc_side == PC_SYMMETRIC) {
 22:     SETERRQ(2,"No symmetric preconditioning for KSPSYMMLQ");
 23:   }
 24:   KSPDefaultGetWork(ksp,9);
 25:   return(0);
 26: }

 28: int  KSPSolve_SYMMLQ(KSP ksp,int *its)
 29: {
 30:   int          ierr,i,maxit;
 31:   PetscScalar  alpha,malpha,beta,mbeta,ibeta,betaold,beta1,ceta,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
 32:   PetscScalar  c=1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,ms,rho0,rho1,rho2,rho3;
 33:   PetscScalar  mone = -1.0,zero = 0.0,dp = 0.0;
 34:   PetscReal    np,s_prod;
 35:   Vec          X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
 36:   Mat          Amat,Pmat;
 37:   MatStructure pflag;
 38:   KSP_SYMMLQ   *symmlq = (KSP_SYMMLQ*)ksp->data;
 39:   PetscTruth   diagonalscale;

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

 45:   maxit   = ksp->max_it;
 46:   X       = ksp->vec_sol;
 47:   B       = ksp->vec_rhs;
 48:   R       = ksp->work[0];
 49:   Z       = ksp->work[1];
 50:   U       = ksp->work[2];
 51:   V       = ksp->work[3];
 52:   W       = ksp->work[4];
 53:   UOLD    = ksp->work[5];
 54:   VOLD    = ksp->work[6];
 55:   Wbar    = ksp->work[7];
 56: 
 57:   PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);

 59:   ksp->its = 0;

 61:   VecSet(&zero,UOLD);          /* u_old <- zeros;  */
 62:   VecCopy(UOLD,VOLD);          /* v_old <- u_old;  */
 63:   VecCopy(UOLD,W);             /* w     <- u_old;  */
 64:   VecCopy(UOLD,Wbar);          /* w_bar <- u_old;  */
 65:   if (!ksp->guess_zero) {
 66:     KSP_MatMult(ksp,Amat,X,R); /*     r <- b - A*x */
 67:     VecAYPX(&mone,B,R);
 68:   } else {
 69:     VecCopy(B,R);              /*     r <- b (x is 0) */
 70:   }

 72:   KSP_PCApply(ksp,ksp->B,R,Z); /* z  <- B*r       */
 73:   VecDot(R,Z,&dp);             /* dp = r'*z;      */
 74:   if (PetscAbsScalar(dp) < symmlq->haptol) {
 75:     PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),symmlq->haptol);
 76:     dp = 0.0;
 77:   }

 79: #if !defined(PETSC_USE_COMPLEX)
 80:   if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
 81: #endif
 82:   dp = PetscSqrtScalar(dp);
 83:   beta = dp;                         /*  beta <- sqrt(r'*z)  */
 84:   beta1 = beta;
 85:   s_prod = PetscAbsScalar(beta1);

 87:   VecCopy(R,V);  /* v <- r; */
 88:   VecCopy(Z,U);  /* u <- z; */
 89:   ibeta = 1.0 / beta;
 90:   VecScale(&ibeta,V);     /* v <- ibeta*v; */
 91:   VecScale(&ibeta,U);     /* u <- ibeta*u; */
 92:   VecCopy(U,Wbar);        /* w_bar <- u;   */
 93:   VecNorm(Z,NORM_2,&np);      /*   np <- ||z||        */
 94:   (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP);  /* test for convergence */
 95:   if (ksp->reason) {*its =  0; return(0);}
 96:   KSPLogResidualHistory(ksp,np);
 97:   KSPMonitor(ksp,0,np);            /* call any registered monitor routines */
 98:   ksp->rnorm = np;

100:   for (i=0; i<maxit; i++){
101:     ksp->its = i+1;

103:     /*    Update    */
104:     if (ksp->its > 1){
105:       VecCopy(V,VOLD);  /* v_old <- v; */
106:       VecCopy(U,UOLD);  /* u_old <- u; */
107: 
108:       ibeta = 1.0 / beta;
109:       VecCopy(R,V);
110:       VecScale(&ibeta,V); /* v <- ibeta*r; */
111:       VecCopy(Z,U);
112:       VecScale(&ibeta,U); /* u <- ibeta*z; */

114:       VecCopy(Wbar,W);
115:       VecScale(&c,W);
116:       VecAXPY(&s,U,W);   /* w  <- c*w_bar + s*u;    (w_k) */
117:       ms = -s;
118:       VecScale(&ms,Wbar);
119:       VecAXPY(&c,U,Wbar); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
120:       VecAXPY(&ceta,W,X); /* x <- x + ceta * w;       (xL_k)  */

122:       ceta_oold = ceta_old;
123:       ceta_old  = ceta;
124:     }

126:     /*   Lanczos  */
127:     KSP_MatMult(ksp,Amat,U,R);   /*  r     <- Amat*u; */
128:     VecDot(U,R,&alpha);          /*  alpha <- u'*r;   */
129:     KSP_PCApply(ksp,ksp->B,R,Z); /*      z <- B*r;    */

131:     malpha = - alpha;
132:     VecAXPY(&malpha,V,R);     /*  r <- r - alpha* v;  */
133:     VecAXPY(&malpha,U,Z);     /*  z <- z - alpha* u;  */
134:     mbeta = - beta;
135:     VecAXPY(&mbeta,VOLD,R);   /*  r <- r - beta * v_old; */
136:     VecAXPY(&mbeta,UOLD,Z);   /*  z <- z - beta * u_old; */
137:     betaold = beta;                                /* beta_k                  */
138:     VecDot(R,Z,&dp);          /* dp <- r'*z;             */
139:     if (PetscAbsScalar(dp) < symmlq->haptol) {
140:       PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),symmlq->haptol);
141:       dp = 0.0;
142:     }

144: #if !defined(PETSC_USE_COMPLEX)
145:      if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
146: #endif
147:      beta = PetscSqrtScalar(dp);                    /*  beta = sqrt(dp); */

149:      /*    QR factorization    */
150:      coold = cold; cold = c; soold = sold; sold = s;
151:      rho0 = cold * alpha - coold * sold * betaold;    /* gamma_bar */
152:      rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);   /* gamma     */
153:      rho2 = sold * alpha + coold * cold * betaold;    /* delta     */
154:      rho3 = soold * betaold;                          /* epsilon   */

156:      /* Givens rotation: [c -s; s c] (different from the Reference!) */
157:      c = rho0 / rho1; s = beta / rho1;

159:      if (ksp->its==1){
160:        ceta = beta1/rho1;
161:      } else {
162:        ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
163:      }
164: 
165:      s_prod = s_prod*PetscAbsScalar(s);
166:      if (c == 0.0){
167:        np = s_prod*1.e16;
168:      } else {
169:        np = s_prod/PetscAbsScalar(c);       /* residual norm for xc_k (CGNORM) */
170:      }
171:      ksp->rnorm = np;
172:      KSPLogResidualHistory(ksp,np);
173:      KSPMonitor(ksp,i+1,np);
174:      (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
175:      if (ksp->reason) break;
176:   }

178:   /* move to the CG point: xc_(k+1) */
179:   if (c == 0.0){
180:     ceta_bar = ceta*1.e15;
181:   } else {
182:     ceta_bar = ceta/c;
183:   }
184:   VecAXPY(&ceta_bar,Wbar,X); /* x <- x + ceta_bar*w_bar */

186:   if (i == maxit) {
187:     ksp->reason = KSP_DIVERGED_ITS;
188:   }
189:   *its = ksp->its;

191:   return(0);
192: }

194: EXTERN_C_BEGIN
195: int KSPCreate_SYMMLQ(KSP ksp)
196: {
197:   KSP_SYMMLQ *symmlq;


202:   ksp->pc_side                   = PC_LEFT;

204:   ierr           = PetscNew(KSP_SYMMLQ,&symmlq);
205:   symmlq->haptol = 1.e-18;
206:   ksp->data      = (void*)symmlq;

208:   /*
209:        Sets the functions that are associated with this data structure 
210:        (in C++ this is the same as defining virtual functions)
211:   */
212:   ksp->ops->setup                = KSPSetUp_SYMMLQ;
213:   ksp->ops->solve                = KSPSolve_SYMMLQ;
214:   ksp->ops->destroy              = KSPDefaultDestroy;
215:   ksp->ops->setfromoptions       = 0;
216:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
217:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;

219:   return(0);
220: }
221: EXTERN_C_END