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