Actual source code: symmlq.c

  1: #define PETSCKSP_DLL

 3:  #include src/ksp/ksp/kspimpl.h

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

 11: PetscErrorCode KSPSetUp_SYMMLQ(KSP ksp)
 12: {

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

 27: PetscErrorCode  KSPSolve_SYMMLQ(KSP ksp)
 28: {
 30:   PetscInt       i;
 31:   PetscScalar    alpha,malpha,beta,mbeta,ibeta,betaold,beta1,ceta = 0,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:   PCDiagonalScale(ksp->pc,&diagonalscale);
 43:   if (diagonalscale) SETERRQ1(PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",ksp->type_name);

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

 58:   ksp->its = 0;

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

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

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

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

102:   i = 0;
103:   do {
104:     ksp->its = i+1;

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

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

125:       ceta_oold = ceta_old;
126:       ceta_old  = ceta;
127:     }

129:     /*   Lanczos  */
130:     KSP_MatMult(ksp,Amat,U,R);   /*  r     <- Amat*u; */
131:     VecDot(U,R,&alpha);          /*  alpha <- u'*r;   */
132:     KSP_PCApply(ksp,R,Z); /*      z <- B*r;    */

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

147: #if !defined(PETSC_USE_COMPLEX)
148:     if (dp < 0.0) {
149:       ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
150:       break;
151:     }
152: #endif
153:     beta = PetscSqrtScalar(dp);                    /*  beta = sqrt(dp); */

155:     /*    QR factorization    */
156:     coold = cold; cold = c; soold = sold; sold = s;
157:     rho0 = cold * alpha - coold * sold * betaold;    /* gamma_bar */
158:     rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta);   /* gamma     */
159:     rho2 = sold * alpha + coold * cold * betaold;    /* delta     */
160:     rho3 = soold * betaold;                          /* epsilon   */

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

165:     if (ksp->its==1){
166:       ceta = beta1/rho1;
167:     } else {
168:       ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
169:     }
170: 
171:     s_prod = s_prod*PetscAbsScalar(s);
172:     if (c == 0.0){
173:       np = s_prod*1.e16;
174:     } else {
175:       np = s_prod/PetscAbsScalar(c);       /* residual norm for xc_k (CGNORM) */
176:     }
177:     ksp->rnorm = np;
178:     KSPLogResidualHistory(ksp,np);
179:     KSPMonitor(ksp,i+1,np);
180:     (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
181:     if (ksp->reason) break;
182:     i++;
183:   } while (i<ksp->max_it);

185:   /* move to the CG point: xc_(k+1) */
186:   if (c == 0.0){
187:     ceta_bar = ceta*1.e15;
188:   } else {
189:     ceta_bar = ceta/c;
190:   }
191:   VecAXPY(X,ceta_bar,Wbar); /* x <- x + ceta_bar*w_bar */

193:   if (i >= ksp->max_it) {
194:     ksp->reason = KSP_DIVERGED_ITS;
195:   }
196:   return(0);
197: }

199: /*MC
200:      KSPSYMMLQ -  This code implements the SYMMLQ method. 

202:    Options Database Keys:
203: .   see KSPSolve()

205:    Level: beginner

207:    Notes: Reference: Paige & Saunders, 1975.

209: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP
210: M*/
214: PetscErrorCode PETSCKSP_DLLEXPORT KSPCreate_SYMMLQ(KSP ksp)
215: {
216:   KSP_SYMMLQ     *symmlq;

220:   ksp->pc_side                   = PC_LEFT;

222:   PetscNew(KSP_SYMMLQ,&symmlq);
223:   symmlq->haptol = 1.e-18;
224:   ksp->data      = (void*)symmlq;

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