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: #undef __FUNCT__
16: int KSPSetUp_SYMMLQ(KSP ksp)
17: {
21: if (ksp->pc_side == PC_RIGHT) {
22: SETERRQ(2,"No right preconditioning for KSPSYMMLQ");
23: } else if (ksp->pc_side == PC_SYMMETRIC) {
24: SETERRQ(2,"No symmetric preconditioning for KSPSYMMLQ");
25: }
26: KSPDefaultGetWork(ksp,9);
27: return(0);
28: }
30: #undef __FUNCT__
32: int KSPSolve_SYMMLQ(KSP ksp,int *its)
33: {
34: int ierr,i,maxit;
35: PetscScalar alpha,malpha,beta,mbeta,ibeta,betaold,beta1,ceta,ceta_oold = 0.0, ceta_old = 0.0,ceta_bar;
36: PetscScalar c=1.0,cold=1.0,s=0.0,sold=0.0,coold,soold,ms,rho0,rho1,rho2,rho3;
37: PetscScalar mone = -1.0,zero = 0.0,dp = 0.0;
38: PetscReal np,s_prod;
39: Vec X,B,R,Z,U,V,W,UOLD,VOLD,Wbar;
40: Mat Amat,Pmat;
41: MatStructure pflag;
42: KSP_SYMMLQ *symmlq = (KSP_SYMMLQ*)ksp->data;
43: PetscTruth diagonalscale;
46: ierr = PCDiagonalScale(ksp->B,&diagonalscale);
47: if (diagonalscale) SETERRQ1(1,"Krylov method %s does not support diagonal scaling",ksp->type_name);
49: maxit = ksp->max_it;
50: X = ksp->vec_sol;
51: B = ksp->vec_rhs;
52: R = ksp->work[0];
53: Z = ksp->work[1];
54: U = ksp->work[2];
55: V = ksp->work[3];
56: W = ksp->work[4];
57: UOLD = ksp->work[5];
58: VOLD = ksp->work[6];
59: Wbar = ksp->work[7];
60:
61: PCGetOperators(ksp->B,&Amat,&Pmat,&pflag);
63: ksp->its = 0;
65: VecSet(&zero,UOLD); /* u_old <- zeros; */
66: VecCopy(UOLD,VOLD); /* v_old <- u_old; */
67: VecCopy(UOLD,W); /* w <- u_old; */
68: VecCopy(UOLD,Wbar); /* w_bar <- u_old; */
69: if (!ksp->guess_zero) {
70: KSP_MatMult(ksp,Amat,X,R); /* r <- b - A*x */
71: VecAYPX(&mone,B,R);
72: } else {
73: VecCopy(B,R); /* r <- b (x is 0) */
74: }
76: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- B*r */
77: VecDot(R,Z,&dp); /* dp = r'*z; */
78: if (PetscAbsScalar(dp) < symmlq->haptol) {
79: PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),symmlq->haptol);
80: dp = 0.0;
81: }
83: #if !defined(PETSC_USE_COMPLEX)
84: if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
85: #endif
86: dp = PetscSqrtScalar(dp);
87: beta = dp; /* beta <- sqrt(r'*z) */
88: beta1 = beta;
89: s_prod = PetscAbsScalar(beta1);
91: VecCopy(R,V); /* v <- r; */
92: VecCopy(Z,U); /* u <- z; */
93: ibeta = 1.0 / beta;
94: VecScale(&ibeta,V); /* v <- ibeta*v; */
95: VecScale(&ibeta,U); /* u <- ibeta*u; */
96: VecCopy(U,Wbar); /* w_bar <- u; */
97: VecNorm(Z,NORM_2,&np); /* np <- ||z|| */
98: (*ksp->converged)(ksp,0,np,&ksp->reason,ksp->cnvP); /* test for convergence */
99: if (ksp->reason) {*its = 0; return(0);}
100: KSPLogResidualHistory(ksp,np);
101: KSPMonitor(ksp,0,np); /* call any registered monitor routines */
102: ksp->rnorm = np;
104: for (i=0; i<maxit; i++){
105: ksp->its = i+1;
107: /* Update */
108: if (ksp->its > 1){
109: VecCopy(V,VOLD); /* v_old <- v; */
110: VecCopy(U,UOLD); /* u_old <- u; */
111:
112: ibeta = 1.0 / beta;
113: VecCopy(R,V);
114: VecScale(&ibeta,V); /* v <- ibeta*r; */
115: VecCopy(Z,U);
116: VecScale(&ibeta,U); /* u <- ibeta*z; */
118: VecCopy(Wbar,W);
119: VecScale(&c,W);
120: VecAXPY(&s,U,W); /* w <- c*w_bar + s*u; (w_k) */
121: ms = -s;
122: VecScale(&ms,Wbar);
123: VecAXPY(&c,U,Wbar); /* w_bar <- -s*w_bar + c*u; (w_bar_(k+1)) */
124: VecAXPY(&ceta,W,X); /* x <- x + ceta * w; (xL_k) */
126: ceta_oold = ceta_old;
127: ceta_old = ceta;
128: }
130: /* Lanczos */
131: KSP_MatMult(ksp,Amat,U,R); /* r <- Amat*u; */
132: VecDot(U,R,&alpha); /* alpha <- u'*r; */
133: KSP_PCApply(ksp,ksp->B,R,Z); /* z <- B*r; */
135: malpha = - alpha;
136: VecAXPY(&malpha,V,R); /* r <- r - alpha* v; */
137: VecAXPY(&malpha,U,Z); /* z <- z - alpha* u; */
138: mbeta = - beta;
139: VecAXPY(&mbeta,VOLD,R); /* r <- r - beta * v_old; */
140: VecAXPY(&mbeta,UOLD,Z); /* z <- z - beta * u_old; */
141: betaold = beta; /* beta_k */
142: VecDot(R,Z,&dp); /* dp <- r'*z; */
143: if (PetscAbsScalar(dp) < symmlq->haptol) {
144: PetscLogInfo(ksp,"KSPSolve_SYMMLQ:Detected happy breakdown %g tolerance %gn",PetscAbsScalar(dp),symmlq->haptol);
145: dp = 0.0;
146: }
148: #if !defined(PETSC_USE_COMPLEX)
149: if (dp < 0.0) SETERRQ(PETSC_ERR_KSP_BRKDWN,"Indefinite preconditioner");
150: #endif
151: beta = PetscSqrtScalar(dp); /* beta = sqrt(dp); */
153: /* QR factorization */
154: coold = cold; cold = c; soold = sold; sold = s;
155: rho0 = cold * alpha - coold * sold * betaold; /* gamma_bar */
156: rho1 = PetscSqrtScalar(rho0*rho0 + beta*beta); /* gamma */
157: rho2 = sold * alpha + coold * cold * betaold; /* delta */
158: rho3 = soold * betaold; /* epsilon */
160: /* Givens rotation: [c -s; s c] (different from the Reference!) */
161: c = rho0 / rho1; s = beta / rho1;
163: if (ksp->its==1){
164: ceta = beta1/rho1;
165: } else {
166: ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
167: }
168:
169: s_prod = s_prod*PetscAbsScalar(s);
170: if (c == 0.0){
171: np = s_prod*1.e16;
172: } else {
173: np = s_prod/PetscAbsScalar(c); /* residual norm for xc_k (CGNORM) */
174: }
175: ksp->rnorm = np;
176: KSPLogResidualHistory(ksp,np);
177: KSPMonitor(ksp,i+1,np);
178: (*ksp->converged)(ksp,i+1,np,&ksp->reason,ksp->cnvP); /* test for convergence */
179: if (ksp->reason) break;
180: }
182: /* move to the CG point: xc_(k+1) */
183: if (c == 0.0){
184: ceta_bar = ceta*1.e15;
185: } else {
186: ceta_bar = ceta/c;
187: }
188: VecAXPY(&ceta_bar,Wbar,X); /* x <- x + ceta_bar*w_bar */
190: if (i == maxit) {
191: ksp->reason = KSP_DIVERGED_ITS;
192: }
193: *its = ksp->its;
195: return(0);
196: }
198: EXTERN_C_BEGIN
199: #undef __FUNCT__
201: int KSPCreate_SYMMLQ(KSP ksp)
202: {
203: KSP_SYMMLQ *symmlq;
208: ksp->pc_side = PC_LEFT;
210: ierr = PetscNew(KSP_SYMMLQ,&symmlq);
211: symmlq->haptol = 1.e-18;
212: ksp->data = (void*)symmlq;
214: /*
215: Sets the functions that are associated with this data structure
216: (in C++ this is the same as defining virtual functions)
217: */
218: ksp->ops->setup = KSPSetUp_SYMMLQ;
219: ksp->ops->solve = KSPSolve_SYMMLQ;
220: ksp->ops->destroy = KSPDefaultDestroy;
221: ksp->ops->setfromoptions = 0;
222: ksp->ops->buildsolution = KSPDefaultBuildSolution;
223: ksp->ops->buildresidual = KSPDefaultBuildResidual;
225: return(0);
226: }
227: EXTERN_C_END