Actual source code: tfqmr.c
petsc-3.3-p7 2013-05-11
2: #include <petsc-private/kspimpl.h>
6: static PetscErrorCode KSPSetUp_TFQMR(KSP ksp)
7: {
11: if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"no symmetric preconditioning for KSPTFQMR");
12: KSPDefaultGetWork(ksp,9);
13: return(0);
14: }
18: static PetscErrorCode KSPSolve_TFQMR(KSP ksp)
19: {
21: PetscInt i,m;
22: PetscScalar rho,rhoold,a,s,b,eta,etaold,psiold,cf;
23: PetscReal dp,dpold,w,dpest,tau,psi,cm;
24: Vec X,B,V,P,R,RP,T,T1,Q,U,D,AUQ;
27: X = ksp->vec_sol;
28: B = ksp->vec_rhs;
29: R = ksp->work[0];
30: RP = ksp->work[1];
31: V = ksp->work[2];
32: T = ksp->work[3];
33: Q = ksp->work[4];
34: P = ksp->work[5];
35: U = ksp->work[6];
36: D = ksp->work[7];
37: T1 = ksp->work[8];
38: AUQ = V;
40: /* Compute initial preconditioned residual */
41: KSPInitialResidual(ksp,X,V,T,R,B);
43: /* Test for nothing to do */
44: VecNorm(R,NORM_2,&dp);
45: PetscObjectTakeAccess(ksp);
46: ksp->rnorm = dp;
47: ksp->its = 0;
48: PetscObjectGrantAccess(ksp);
49: KSPMonitor(ksp,0,dp);
50: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
51: if (ksp->reason) return(0);
53: /* Make the initial Rp == R */
54: VecCopy(R,RP);
56: /* Set the initial conditions */
57: etaold = 0.0;
58: psiold = 0.0;
59: tau = dp;
60: dpold = dp;
62: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
63: VecCopy(R,U);
64: VecCopy(R,P);
65: KSP_PCApplyBAorAB(ksp,P,V,T);
66: VecSet(D,0.0);
68: i=0;
69: do {
70: PetscObjectTakeAccess(ksp);
71: ksp->its++;
72: PetscObjectGrantAccess(ksp);
73: VecDot(V,RP,&s); /* s <- (v,rp) */
74: a = rhoold / s; /* a <- rho / s */
75: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
76: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
77: KSP_PCApplyBAorAB(ksp,T,AUQ,T1);
78: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
79: VecNorm(R,NORM_2,&dp);
80: for (m=0; m<2; m++) {
81: if (!m) {
82: w = PetscSqrtReal(dp*dpold);
83: } else {
84: w = dp;
85: }
86: psi = w / tau;
87: cm = 1.0 / PetscSqrtReal(1.0 + psi * psi);
88: tau = tau * psi * cm;
89: eta = cm * cm * a;
90: cf = psiold * psiold * etaold / a;
91: if (!m) {
92: VecAYPX(D,cf,U);
93: } else {
94: VecAYPX(D,cf,Q);
95: }
96: VecAXPY(X,eta,D);
98: dpest = PetscSqrtReal(m + 1.0) * tau;
99: PetscObjectTakeAccess(ksp);
100: ksp->rnorm = dpest;
101: PetscObjectGrantAccess(ksp);
102: KSPLogResidualHistory(ksp,dpest);
103: KSPMonitor(ksp,i+1,dpest);
104: (*ksp->converged)(ksp,i+1,dpest,&ksp->reason,ksp->cnvP);
105: if (ksp->reason) break;
107: etaold = eta;
108: psiold = psi;
109: }
110: if (ksp->reason) break;
112: VecDot(R,RP,&rho); /* rho <- (r,rp) */
113: b = rho / rhoold; /* b <- rho / rhoold */
114: VecWAXPY(U,b,Q,R); /* u <- r + b q */
115: VecAXPY(Q,b,P);
116: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
117: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
119: rhoold = rho;
120: dpold = dp;
122: i++;
123: } while (i<ksp->max_it);
124: if (i >= ksp->max_it) {
125: ksp->reason = KSP_DIVERGED_ITS;
126: }
128: KSPUnwindPreconditioner(ksp,X,T);
129: return(0);
130: }
132: /*MC
133: KSPTFQMR - A transpose free QMR (quasi minimal residual),
135: Options Database Keys:
136: . see KSPSolve()
138: Level: beginner
140: Notes: Supports left and right preconditioning, but not symmetric
142: The "residual norm" computed in this algorithm is actually just an upper bound on the actual residual norm.
143: That is for left preconditioning it is a bound on the preconditioned residual and for right preconditioning
144: it is a bound on the true residual.
146: References: Freund, 1993
148: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTCQMR
149: M*/
150: EXTERN_C_BEGIN
153: PetscErrorCode KSPCreate_TFQMR(KSP ksp)
154: {
158: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
159: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);
160: ksp->data = (void*)0;
161: ksp->ops->setup = KSPSetUp_TFQMR;
162: ksp->ops->solve = KSPSolve_TFQMR;
163: ksp->ops->destroy = KSPDefaultDestroy;
164: ksp->ops->buildsolution = KSPDefaultBuildSolution;
165: ksp->ops->buildresidual = KSPDefaultBuildResidual;
166: ksp->ops->setfromoptions = 0;
167: ksp->ops->view = 0;
168: return(0);
169: }
170: EXTERN_C_END