Actual source code: tcqmr.c
petsc-3.9.2 2018-05-20
2: /*
3: This file contains an implementation of Tony Chan's transpose-free QMR.
5: Note: The vector dot products in the code have not been checked for the
6: complex numbers version, so most probably some are incorrect.
7: */
9: #include <../src/ksp/ksp/impls/tcqmr/tcqmrimpl.h>
11: static PetscErrorCode KSPSolve_TCQMR(KSP ksp)
12: {
13: PetscReal rnorm0,rnorm,dp1,Gamma;
14: PetscScalar theta,ep,cl1,sl1,cl,sl,sprod,tau_n1,f;
15: PetscScalar deltmp,rho,beta,eptmp,ta,s,c,tau_n,delta;
16: PetscScalar dp11,dp2,rhom1,alpha,tmp;
20: ksp->its = 0;
22: KSPInitialResidual(ksp,x,u,v,r,b);
23: VecNorm(r,NORM_2,&rnorm0); /* rnorm0 = ||r|| */
25: (*ksp->converged)(ksp,0,rnorm0,&ksp->reason,ksp->cnvP);
26: if (ksp->reason) return(0);
28: VecSet(um1,0.0);
29: VecCopy(r,u);
30: rnorm = rnorm0;
31: tmp = 1.0/rnorm; VecScale(u,tmp);
32: VecSet(vm1,0.0);
33: VecCopy(u,v);
34: VecCopy(u,v0);
35: VecSet(pvec1,0.0);
36: VecSet(pvec2,0.0);
37: VecSet(p,0.0);
38: theta = 0.0;
39: ep = 0.0;
40: cl1 = 0.0;
41: sl1 = 0.0;
42: cl = 0.0;
43: sl = 0.0;
44: sprod = 1.0;
45: tau_n1= rnorm0;
46: f = 1.0;
47: Gamma = 1.0;
48: rhom1 = 1.0;
50: /*
51: CALCULATE SQUARED LANCZOS vectors
52: */
53: (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
54: while (!ksp->reason) {
55: KSPMonitor(ksp,ksp->its,rnorm);
56: ksp->its++;
58: KSP_PCApplyBAorAB(ksp,u,y,vtmp); /* y = A*u */
59: VecDot(y,v0,&dp11);
60: VecDot(u,v0,&dp2);
61: alpha = dp11 / dp2; /* alpha = v0'*y/v0'*u */
62: deltmp = alpha;
63: VecCopy(y,z);
64: VecAXPY(z,-alpha,u); /* z = y - alpha u */
65: VecDot(u,v0,&rho);
66: beta = rho / (f*rhom1);
67: rhom1 = rho;
68: VecCopy(z,utmp); /* up1 = (A-alpha*I)*
69: (z-2*beta*p) + f*beta*
70: beta*um1 */
71: VecAXPY(utmp,-2.0*beta,p);
72: KSP_PCApplyBAorAB(ksp,utmp,up1,vtmp);
73: VecAXPY(up1,-alpha,utmp);
74: VecAXPY(up1,f*beta*beta,um1);
75: VecNorm(up1,NORM_2,&dp1);
76: f = 1.0 / dp1;
77: VecScale(up1,f);
78: VecAYPX(p,-beta,z); /* p = f*(z-beta*p) */
79: VecScale(p,f);
80: VecCopy(u,um1);
81: VecCopy(up1,u);
82: beta = beta/Gamma;
83: eptmp = beta;
84: KSP_PCApplyBAorAB(ksp,v,vp1,vtmp);
85: VecAXPY(vp1,-alpha,v);
86: VecAXPY(vp1,-beta,vm1);
87: VecNorm(vp1,NORM_2,&Gamma);
88: VecScale(vp1,1.0/Gamma);
89: VecCopy(v,vm1);
90: VecCopy(vp1,v);
92: /*
93: SOLVE Ax = b
94: */
95: /* Apply last two Given's (Gl-1 and Gl) rotations to (beta,alpha,Gamma) */
96: if (ksp->its > 2) {
97: theta = sl1*beta;
98: eptmp = -cl1*beta;
99: }
100: if (ksp->its > 1) {
101: ep = -cl*eptmp + sl*alpha;
102: deltmp = -sl*eptmp - cl*alpha;
103: }
104: if (PetscAbsReal(Gamma) > PetscAbsScalar(deltmp)) {
105: ta = -deltmp / Gamma;
106: s = 1.0 / PetscSqrtScalar(1.0 + ta*ta);
107: c = s*ta;
108: } else {
109: ta = -Gamma/deltmp;
110: c = 1.0 / PetscSqrtScalar(1.0 + ta*ta);
111: s = c*ta;
112: }
114: delta = -c*deltmp + s*Gamma;
115: tau_n = -c*tau_n1; tau_n1 = -s*tau_n1;
116: VecCopy(vm1,pvec);
117: VecAXPY(pvec,-theta,pvec2);
118: VecAXPY(pvec,-ep,pvec1);
119: VecScale(pvec,1.0/delta);
120: VecAXPY(x,tau_n,pvec);
121: cl1 = cl; sl1 = sl; cl = c; sl = s;
123: VecCopy(pvec1,pvec2);
124: VecCopy(pvec,pvec1);
126: /* Compute the upper bound on the residual norm r (See QMR paper p. 13) */
127: sprod = sprod*PetscAbsScalar(s);
128: rnorm = rnorm0 * PetscSqrtReal((PetscReal)ksp->its+2.0) * PetscRealPart(sprod);
129: (*ksp->converged)(ksp,ksp->its,rnorm,&ksp->reason,ksp->cnvP);
130: if (ksp->its >= ksp->max_it) {
131: if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS;
132: break;
133: }
134: }
135: KSPMonitor(ksp,ksp->its,rnorm);
136: KSPUnwindPreconditioner(ksp,x,vtmp);
137: return(0);
138: }
140: static PetscErrorCode KSPSetUp_TCQMR(KSP ksp)
141: {
145: if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"no symmetric preconditioning for KSPTCQMR");
146: KSPSetWorkVecs(ksp,TCQMR_VECS);
147: return(0);
148: }
150: /*MC
151: KSPTCQMR - A variant of QMR (quasi minimal residual) developed by Tony Chan
153: Options Database Keys:
154: . see KSPSolve()
156: Level: beginner
158: Notes: Supports either left or right preconditioning, but not symmetric
160: The "residual norm" computed in this algorithm is actually just an upper bound on the actual residual norm.
161: That is for left preconditioning it is a bound on the preconditioned residual and for right preconditioning
162: it is a bound on the true residual.
164: References:
165: . 1. - Tony F. Chan, Lisette de Pillis, and Henk van der Vorst, Transpose free formulations of Lanczos type methods for nonsymmetric linear systems,
166: Numerical Algorithms, Volume 17, 1998.
168: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTFQMR
170: M*/
172: PETSC_EXTERN PetscErrorCode KSPCreate_TCQMR(KSP ksp)
173: {
177: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
178: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
180: ksp->data = (void*)0;
181: ksp->ops->buildsolution = KSPBuildSolutionDefault;
182: ksp->ops->buildresidual = KSPBuildResidualDefault;
183: ksp->ops->setup = KSPSetUp_TCQMR;
184: ksp->ops->solve = KSPSolve_TCQMR;
185: ksp->ops->destroy = KSPDestroyDefault;
186: ksp->ops->setfromoptions = 0;
187: ksp->ops->view = 0;
188: return(0);
189: }