Actual source code: tfqmr.c
1: /*$Id: tfqmr.c,v 1.61 2001/08/07 03:03:54 balay Exp $*/
3: /*
4: This code implements the TFQMR (Transpose-free variant of Quasi-Minimal
5: Residual) method. Reference: Freund, 1993
7: Note that for the complex numbers version, the VecDot() arguments
8: within the code MUST remain in the order given for correct computation
9: of inner products.
10: */
12: #include src/sles/ksp/kspimpl.h
14: static int KSPSetUp_TFQMR(KSP ksp)
15: {
19: if (ksp->pc_side == PC_SYMMETRIC){
20: SETERRQ(2,"no symmetric preconditioning for KSPTFQMR");
21: }
22: KSPDefaultGetWork(ksp,9);
23: return(0);
24: }
26: static int KSPSolve_TFQMR(KSP ksp,int *its)
27: {
28: int i,maxit,m, ierr;
29: PetscScalar rho,rhoold,a,s,b,eta,etaold,psiold,cf,tmp,one = 1.0,zero = 0.0;
30: PetscReal dp,dpold,w,dpest,tau,psi,cm;
31: Vec X,B,V,P,R,RP,T,T1,Q,U,D,AUQ;
34: maxit = ksp->max_it;
35: X = ksp->vec_sol;
36: B = ksp->vec_rhs;
37: R = ksp->work[0];
38: RP = ksp->work[1];
39: V = ksp->work[2];
40: T = ksp->work[3];
41: Q = ksp->work[4];
42: P = ksp->work[5];
43: U = ksp->work[6];
44: D = ksp->work[7];
45: T1 = ksp->work[8];
46: AUQ = V;
48: /* Compute initial preconditioned residual */
49: KSPInitialResidual(ksp,X,V,T,R,B);
51: /* Test for nothing to do */
52: VecNorm(R,NORM_2,&dp);
53: PetscObjectTakeAccess(ksp);
54: ksp->rnorm = dp;
55: ksp->its = 0;
56: PetscObjectGrantAccess(ksp);
57: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
58: if (ksp->reason) {*its = 0; return(0);}
59: KSPMonitor(ksp,0,dp);
61: /* Make the initial Rp == R */
62: VecCopy(R,RP);
64: /* Set the initial conditions */
65: etaold = 0.0;
66: psiold = 0.0;
67: tau = dp;
68: dpold = dp;
70: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
71: VecCopy(R,U);
72: VecCopy(R,P);
73: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T);
74: VecSet(&zero,D);
76: for (i=0; i<maxit; i++) {
77: PetscObjectTakeAccess(ksp);
78: ksp->its++;
79: PetscObjectGrantAccess(ksp);
80: VecDot(V,RP,&s); /* s <- (v,rp) */
81: a = rhoold / s; /* a <- rho / s */
82: tmp = -a; VecWAXPY(&tmp,V,U,Q); /* q <- u - a v */
83: VecWAXPY(&one,U,Q,T); /* t <- u + q */
84: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,T,AUQ,T1);
85: VecAXPY(&tmp,AUQ,R); /* r <- r - a K (u + q) */
86: VecNorm(R,NORM_2,&dp);
87: for (m=0; m<2; m++) {
88: if (!m) {
89: w = sqrt(dp*dpold);
90: } else {
91: w = dp;
92: }
93: psi = w / tau;
94: cm = 1.0 / sqrt(1.0 + psi * psi);
95: tau = tau * psi * cm;
96: eta = cm * cm * a;
97: cf = psiold * psiold * etaold / a;
98: if (!m) {
99: VecAYPX(&cf,U,D);
100: } else {
101: VecAYPX(&cf,Q,D);
102: }
103: VecAXPY(&eta,D,X);
105: dpest = sqrt(m + 1.0) * tau;
106: PetscObjectTakeAccess(ksp);
107: ksp->rnorm = dpest;
108: PetscObjectGrantAccess(ksp);
109: KSPLogResidualHistory(ksp,dpest);
110: KSPMonitor(ksp,i+1,dpest);
111: (*ksp->converged)(ksp,i+1,dpest,&ksp->reason,ksp->cnvP);
112: if (ksp->reason) break;
114: etaold = eta;
115: psiold = psi;
116: }
117: if (ksp->reason) break;
119: VecDot(R,RP,&rho); /* rho <- (r,rp) */
120: b = rho / rhoold; /* b <- rho / rhoold */
121: VecWAXPY(&b,Q,R,U); /* u <- r + b q */
122: VecAXPY(&b,P,Q);
123: VecWAXPY(&b,Q,U,P); /* p <- u + b(q + b p) */
124: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,Q); /* v <- K p */
126: rhoold = rho;
127: dpold = dp;
128: }
129: if (i == maxit) {
130: i--;
131: ksp->reason = KSP_DIVERGED_ITS;
132: }
134: KSPUnwindPreconditioner(ksp,X,T);
135: *its = i + 1;
136: return(0);
137: }
139: EXTERN_C_BEGIN
140: int KSPCreate_TFQMR(KSP ksp)
141: {
143: ksp->data = (void*)0;
144: ksp->pc_side = PC_LEFT;
145: ksp->ops->setup = KSPSetUp_TFQMR;
146: ksp->ops->solve = KSPSolve_TFQMR;
147: ksp->ops->destroy = KSPDefaultDestroy;
148: ksp->ops->buildsolution = KSPDefaultBuildSolution;
149: ksp->ops->buildresidual = KSPDefaultBuildResidual;
150: ksp->ops->setfromoptions = 0;
151: ksp->ops->view = 0;
152: return(0);
153: }
154: EXTERN_C_END