Actual source code: symmlq.c

petsc-3.3-p7 2013-05-11
  2: #include <petsc-private/kspimpl.h>

  4: typedef struct {
  5:   PetscReal haptol;
  6: } KSP_SYMMLQ;

 10: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
 11: {

 15:   KSPDefaultGetWork(ksp,9);
 16:   return(0);
 17: }

 21: PetscErrorCode  KSPSolve_SYMMLQ(KSP ksp)
 22: {
 24:   PetscInt       i;
 25:   PetscScalar    alpha,beta,ibeta,betaold,beta1,ceta = 0,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
 26:   PetscScalar    c=1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,rho0,rho1,rho2,rho3;
 27:   PetscScalar    dp = 0.0;
 28:   PetscReal      np,s_prod;
 29:   Vec            X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
 30:   Mat            Amat,Pmat;
 31:   MatStructure   pflag;
 32:   KSP_SYMMLQ     *symmlq = (KSP_SYMMLQ*)ksp->data;
 33:   PetscBool      diagonalscale;

 36:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
 37:   if (diagonalscale) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

 39:   X       = ksp->vec_sol;
 40:   B       = ksp->vec_rhs;
 41:   R       = ksp->work[0];
 42:   Z       = ksp->work[1];
 43:   U       = ksp->work[2];
 44:   V       = ksp->work[3];
 45:   W       = ksp->work[4];
 46:   UOLD    = ksp->work[5];
 47:   VOLD    = ksp->work[6];
 48:   Wbar    = ksp->work[7];
 49: 
 50:   PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);

 52:   ksp->its = 0;

 54:   VecSet(UOLD,0.0);           /* u_old <- zeros;  */
 55:   VecCopy(UOLD,VOLD);          /* v_old <- u_old;  */
 56:   VecCopy(UOLD,W);             /* w     <- u_old;  */
 57:   VecCopy(UOLD,Wbar);          /* w_bar <- u_old;  */
 58:   if (!ksp->guess_zero) {
 59:     KSP_MatMult(ksp,Amat,X,R); /*     r <- b - A*x */
 60:     VecAYPX(R,-1.0,B);
 61:   } else {
 62:     VecCopy(B,R);              /*     r <- b (x is 0) */
 63:   }

 65:   KSP_PCApply(ksp,R,Z); /* z  <- B*r       */
 66:   VecDot(R,Z,&dp);             /* dp = r'*z;      */
 67:   if (PetscAbsScalar(dp) < symmlq->haptol) {
 68:     PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),symmlq->haptol);
 69:     ksp->rnorm  = 0.0;  /* what should we really put here? */
 70:     ksp->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;  /* bugfix proposed by Lourens (lourens.vanzanen@shell.com) */
 71:     return(0);
 72:   }

 74: #if !defined(PETSC_USE_COMPLEX)
 75:   if (dp < 0.0) {
 76:     ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
 77:     return(0);
 78:   }
 79: #endif
 80:   dp     = PetscSqrtScalar(dp);
 81:   beta   = dp;                         /*  beta <- sqrt(r'*z)  */
 82:   beta1  = beta;
 83:   s_prod = PetscAbsScalar(beta1);

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

 98:   i = 0; ceta = 0.;
 99:   do {
100:     ksp->its = i+1;

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

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

119:       ceta_oold = ceta_old;
120:       ceta_old  = ceta;
121:     }

123:     /*   Lanczos  */
124:     KSP_MatMult(ksp,Amat,U,R);   /*  r     <- Amat*u; */
125:     VecDot(U,R,&alpha);          /*  alpha <- u'*r;   */
126:     KSP_PCApply(ksp,R,Z); /*      z <- B*r;    */

128:     VecAXPY(R,-alpha,V);     /*  r <- r - alpha* v;  */
129:     VecAXPY(Z,-alpha,U);     /*  z <- z - alpha* u;  */
130:     VecAXPY(R,-beta,VOLD);   /*  r <- r - beta * v_old; */
131:     VecAXPY(Z,-beta,UOLD);   /*  z <- z - beta * u_old; */
132:     betaold = beta;                                /* beta_k                  */
133:     VecDot(R,Z,&dp);          /* dp <- r'*z;             */
134:     if (PetscAbsScalar(dp) < symmlq->haptol) {
135:       PetscInfo2(ksp,"Detected happy breakdown %G tolerance %G\n",PetscAbsScalar(dp),symmlq->haptol);
136:       dp = 0.0;
137:     }

139: #if !defined(PETSC_USE_COMPLEX)
140:     if (dp < 0.0) {
141:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
142:       break;
143:     }
144: #endif
145:     beta = PetscSqrtScalar(dp);                    /*  beta = sqrt(dp); */

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

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

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

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

185:   if (i >= ksp->max_it) {
186:     ksp->reason = KSP_DIVERGED_ITS;
187:   }
188:   return(0);
189: }

191: /*MC
192:      KSPSYMMLQ -  This code implements the SYMMLQ method. 

194:    Options Database Keys:
195: .   see KSPSolve()

197:    Level: beginner

199:    Notes: The operator and the preconditioner must be symmetric for this method. The 
200:           preconditioner must be POSITIVE-DEFINITE.

202:           Supports only left preconditioning.

204:    Reference: Paige & Saunders, 1975.

206: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
207: M*/
208: EXTERN_C_BEGIN
211: PetscErrorCode  KSPCreate_SYMMLQ(KSP ksp)
212: {
213:   KSP_SYMMLQ     *symmlq;

217:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
218:   PetscNewLog(ksp,KSP_SYMMLQ,&symmlq);
219:   symmlq->haptol = 1.e-18;
220:   ksp->data      = (void*)symmlq;

222:   /*
223:        Sets the functions that are associated with this data structure 
224:        (in C++ this is the same as defining virtual functions)
225:   */
226:   ksp->ops->setup                = KSPSetUp_SYMMLQ;
227:   ksp->ops->solve                = KSPSolve_SYMMLQ;
228:   ksp->ops->destroy              = KSPDefaultDestroy;
229:   ksp->ops->setfromoptions       = 0;
230:   ksp->ops->buildsolution        = KSPDefaultBuildSolution;
231:   ksp->ops->buildresidual        = KSPDefaultBuildResidual;
232:   return(0);
233: }
234: EXTERN_C_END