Actual source code: tcqmr.c
1: /*$Id: tcqmr.c,v 1.59 2001/08/07 03:03:54 balay Exp $*/
3: /*
4: This file contains an implementation of Tony Chan's transpose-free QMR.
6: Note: The vector dot products in the code have not been checked for the
7: complex numbers version, so most probably some are incorrect.
8: */
10: #include src/sles/ksp/kspimpl.h
11: #include src/sles/ksp/impls/tcqmr/tcqmrp.h
13: static int KSPSolve_TCQMR(KSP ksp,int *its)
14: {
15: PetscReal rnorm0,rnorm,dp1,Gamma;
16: PetscScalar theta,ep,cl1,sl1,cl,sl,sprod,tau_n1,f;
17: PetscScalar deltmp,rho,beta,eptmp,ta,s,c,tau_n,delta;
18: PetscScalar dp11,dp2,rhom1,alpha,tmp,zero = 0.0;
19: int ierr;
22: ksp->its = 0;
24: ierr = KSPInitialResidual(ksp,x,u,v,r,b);
25: ierr = VecNorm(r,NORM_2,&rnorm0); /* rnorm0 = ||r|| */
27: (*ksp->converged)(ksp,0,rnorm0,&ksp->reason,ksp->cnvP);
28: if (ksp->reason) {*its = 0; return(0);}
30: ierr = VecSet(&zero,um1);
31: ierr = VecCopy(r,u);
32: rnorm = rnorm0;
33: tmp = 1.0/rnorm; VecScale(&tmp,u);
34: ierr = VecSet(&zero,vm1);
35: ierr = VecCopy(u,v);
36: ierr = VecCopy(u,v0);
37: ierr = VecSet(&zero,pvec1);
38: ierr = VecSet(&zero,pvec2);
39: ierr = VecSet(&zero,p);
40: theta = 0.0;
41: ep = 0.0;
42: cl1 = 0.0;
43: sl1 = 0.0;
44: cl = 0.0;
45: sl = 0.0;
46: sprod = 1.0;
47: tau_n1= rnorm0;
48: f = 1.0;
49: Gamma = 1.0;
50: rhom1 = 1.0;
52: /*
53: CALCULATE SQUARED LANCZOS vectors
54: */
55: (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
56: while (!ksp->reason){
57: KSPMonitor(ksp,ksp->its,rnorm);
58: ksp->its++;
60: ierr = KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,u,y,vtmp); /* y = A*u */
61: ierr = VecDot(v0,y,&dp11);
62: ierr = VecDot(v0,u,&dp2);
63: alpha = dp11 / dp2; /* alpha = v0'*y/v0'*u */
64: deltmp = alpha;
65: ierr = VecCopy(y,z);
66: tmp = -alpha;
67: ierr = VecAXPY(&tmp,u,z); /* z = y - alpha u */
68: ierr = VecDot(v0,u,&rho);
69: beta = rho / (f*rhom1);
70: rhom1 = rho;
71: ierr = VecCopy(z,utmp); /* up1 = (A-alpha*I)*
72: (z-2*beta*p) + f*beta*
73: beta*um1 */
74: tmp = -2.0*beta;VecAXPY(&tmp,p,utmp);
75: ierr = KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,utmp,up1,vtmp);
76: tmp = -alpha; VecAXPY(&tmp,utmp,up1);
77: tmp = f*beta*beta; VecAXPY(&tmp,um1,up1);
78: ierr = VecNorm(up1,NORM_2,&dp1);
79: f = 1.0 / dp1;
80: ierr = VecScale(&f,up1);
81: tmp = -beta;
82: ierr = VecAYPX(&tmp,z,p); /* p = f*(z-beta*p) */
83: ierr = VecScale(&f,p);
84: ierr = VecCopy(u,um1);
85: ierr = VecCopy(up1,u);
86: beta = beta/Gamma;
87: eptmp = beta;
88: ierr = KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,v,vp1,vtmp);
89: tmp = -alpha; VecAXPY(&tmp,v,vp1);
90: tmp = -beta; VecAXPY(&tmp,vm1,vp1);
91: ierr = VecNorm(vp1,NORM_2,&Gamma);
92: tmp = 1.0/Gamma; VecScale(&tmp,vp1);
93: ierr = VecCopy(v,vm1);
94: ierr = VecCopy(vp1,v);
96: /*
97: SOLVE Ax = b
98: */
99: /* Apply last two Given's (Gl-1 and Gl) rotations to (beta,alpha,Gamma) */
100: if (ksp->its > 1) {
101: theta = sl1*beta;
102: eptmp = -cl1*beta;
103: }
104: if (ksp->its > 0) {
105: ep = -cl*eptmp + sl*alpha;
106: deltmp = -sl*eptmp - cl*alpha;
107: }
108: if (PetscAbsReal(Gamma) > PetscAbsScalar(deltmp)) {
109: ta = -deltmp / Gamma;
110: s = 1.0 / PetscSqrtScalar(1.0 + ta*ta);
111: c = s*ta;
112: } else {
113: ta = -Gamma/deltmp;
114: c = 1.0 / PetscSqrtScalar(1.0 + ta*ta);
115: s = c*ta;
116: }
118: delta = -c*deltmp + s*Gamma;
119: tau_n = -c*tau_n1; tau_n1 = -s*tau_n1;
120: ierr = VecCopy(vm1,pvec);
121: tmp = -theta; VecAXPY(&tmp,pvec2,pvec);
122: tmp = -ep; VecAXPY(&tmp,pvec1,pvec);
123: tmp = 1.0/delta; VecScale(&tmp,pvec);
124: ierr = VecAXPY(&tau_n,pvec,x);
125: cl1 = cl; sl1 = sl; cl = c; sl = s;
127: VecCopy(pvec1,pvec2);
128: VecCopy(pvec,pvec1);
130: /* Compute the upper bound on the residual norm r (See QMR paper p. 13) */
131: sprod = sprod*PetscAbsScalar(s);
132: #if defined(PETSC_USE_COMPLEX)
133: rnorm = rnorm0 * sqrt((double)ksp->its+2.0) * PetscRealPart(sprod);
134: #else
135: rnorm = rnorm0 * sqrt((double)ksp->its+2.0) * sprod;
136: #endif
137: if (ksp->its >= ksp->max_it) {ksp->reason = KSP_DIVERGED_ITS; break;}
138: }
139: KSPMonitor(ksp,ksp->its,rnorm);
140: KSPUnwindPreconditioner(ksp,x,vtmp);
142: *its = ksp->its;
143: return(0);
144: }
146: static int KSPSetUp_TCQMR(KSP ksp)
147: {
151: if (ksp->pc_side == PC_SYMMETRIC){
152: SETERRQ(2,"no symmetric preconditioning for KSPTCQMR");
153: }
154: KSPDefaultGetWork(ksp,TCQMR_VECS);
155: return(0);
156: }
158: EXTERN_C_BEGIN
159: int KSPCreate_TCQMR(KSP ksp)
160: {
162: ksp->data = (void*)0;
163: ksp->pc_side = PC_LEFT;
164: ksp->ops->buildsolution = KSPDefaultBuildSolution;
165: ksp->ops->buildresidual = KSPDefaultBuildResidual;
166: ksp->ops->setup = KSPSetUp_TCQMR;
167: ksp->ops->solve = KSPSolve_TCQMR;
168: ksp->ops->destroy = KSPDefaultDestroy;
169: ksp->ops->setfromoptions = 0;
170: ksp->ops->view = 0;
171: return(0);
172: }
173: EXTERN_C_END