Actual source code: tfqmr.c
1: /*$Id: tfqmr.c,v 1.61 2001/08/07 03:03:54 balay Exp $*/
3: #include src/ksp/ksp/kspimpl.h
7: static int KSPSetUp_TFQMR(KSP ksp)
8: {
12: if (ksp->pc_side == PC_SYMMETRIC){
13: SETERRQ(2,"no symmetric preconditioning for KSPTFQMR");
14: }
15: KSPDefaultGetWork(ksp,9);
16: return(0);
17: }
21: static int KSPSolve_TFQMR(KSP ksp)
22: {
23: int i,m, ierr;
24: PetscScalar rho,rhoold,a,s,b,eta,etaold,psiold,cf,tmp,one = 1.0,zero = 0.0;
25: PetscReal dp,dpold,w,dpest,tau,psi,cm;
26: Vec X,B,V,P,R,RP,T,T1,Q,U,D,AUQ;
29: X = ksp->vec_sol;
30: B = ksp->vec_rhs;
31: R = ksp->work[0];
32: RP = ksp->work[1];
33: V = ksp->work[2];
34: T = ksp->work[3];
35: Q = ksp->work[4];
36: P = ksp->work[5];
37: U = ksp->work[6];
38: D = ksp->work[7];
39: T1 = ksp->work[8];
40: AUQ = V;
42: /* Compute initial preconditioned residual */
43: KSPInitialResidual(ksp,X,V,T,R,B);
45: /* Test for nothing to do */
46: VecNorm(R,NORM_2,&dp);
47: PetscObjectTakeAccess(ksp);
48: ksp->rnorm = dp;
49: ksp->its = 0;
50: PetscObjectGrantAccess(ksp);
51: KSPMonitor(ksp,0,dp);
52: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
53: if (ksp->reason) return(0);
55: /* Make the initial Rp == R */
56: VecCopy(R,RP);
58: /* Set the initial conditions */
59: etaold = 0.0;
60: psiold = 0.0;
61: tau = dp;
62: dpold = dp;
64: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
65: VecCopy(R,U);
66: VecCopy(R,P);
67: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,T);
68: VecSet(&zero,D);
70: i=0;
71: do {
72: PetscObjectTakeAccess(ksp);
73: ksp->its++;
74: PetscObjectGrantAccess(ksp);
75: VecDot(V,RP,&s); /* s <- (v,rp) */
76: a = rhoold / s; /* a <- rho / s */
77: tmp = -a; VecWAXPY(&tmp,V,U,Q); /* q <- u - a v */
78: VecWAXPY(&one,U,Q,T); /* t <- u + q */
79: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,T,AUQ,T1);
80: VecAXPY(&tmp,AUQ,R); /* r <- r - a K (u + q) */
81: VecNorm(R,NORM_2,&dp);
82: for (m=0; m<2; m++) {
83: if (!m) {
84: w = sqrt(dp*dpold);
85: } else {
86: w = dp;
87: }
88: psi = w / tau;
89: cm = 1.0 / sqrt(1.0 + psi * psi);
90: tau = tau * psi * cm;
91: eta = cm * cm * a;
92: cf = psiold * psiold * etaold / a;
93: if (!m) {
94: VecAYPX(&cf,U,D);
95: } else {
96: VecAYPX(&cf,Q,D);
97: }
98: VecAXPY(&eta,D,X);
100: dpest = sqrt(m + 1.0) * tau;
101: PetscObjectTakeAccess(ksp);
102: ksp->rnorm = dpest;
103: PetscObjectGrantAccess(ksp);
104: KSPLogResidualHistory(ksp,dpest);
105: KSPMonitor(ksp,i+1,dpest);
106: (*ksp->converged)(ksp,i+1,dpest,&ksp->reason,ksp->cnvP);
107: if (ksp->reason) break;
109: etaold = eta;
110: psiold = psi;
111: }
112: if (ksp->reason) break;
114: VecDot(R,RP,&rho); /* rho <- (r,rp) */
115: b = rho / rhoold; /* b <- rho / rhoold */
116: VecWAXPY(&b,Q,R,U); /* u <- r + b q */
117: VecAXPY(&b,P,Q);
118: VecWAXPY(&b,Q,U,P); /* p <- u + b(q + b p) */
119: KSP_PCApplyBAorAB(ksp,ksp->B,ksp->pc_side,P,V,Q); /* v <- K p */
121: rhoold = rho;
122: dpold = dp;
124: i++;
125: } while (i<ksp->max_it);
126: if (i == ksp->max_it) {
127: ksp->reason = KSP_DIVERGED_ITS;
128: }
130: KSPUnwindPreconditioner(ksp,X,T);
131: return(0);
132: }
134: /*MC
135: KSPRTFQMR - A transpose free QMR (quasi minimal residual), Freund, 1993
137: Options Database Keys:
138: . see KSPSolve()
140: Level: beginner
142: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTCQMR
144: M*/
146: EXTERN_C_BEGIN
149: int KSPCreate_TFQMR(KSP ksp)
150: {
152: ksp->data = (void*)0;
153: ksp->pc_side = PC_LEFT;
154: ksp->ops->setup = KSPSetUp_TFQMR;
155: ksp->ops->solve = KSPSolve_TFQMR;
156: ksp->ops->destroy = KSPDefaultDestroy;
157: ksp->ops->buildsolution = KSPDefaultBuildSolution;
158: ksp->ops->buildresidual = KSPDefaultBuildResidual;
159: ksp->ops->setfromoptions = 0;
160: ksp->ops->view = 0;
161: return(0);
162: }
163: EXTERN_C_END