Actual source code: tfqmr.c
petsc-3.9.1 2018-04-29
2: #include <petsc/private/kspimpl.h>
4: static PetscErrorCode KSPSetUp_TFQMR(KSP ksp)
5: {
9: if (ksp->pc_side == PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"no symmetric preconditioning for KSPTFQMR");
10: KSPSetWorkVecs(ksp,9);
11: return(0);
12: }
14: static PetscErrorCode KSPSolve_TFQMR(KSP ksp)
15: {
17: PetscInt i,m;
18: PetscScalar rho,rhoold,a,s,b,eta,etaold,psiold,cf;
19: PetscReal dp,dpold,w,dpest,tau,psi,cm;
20: Vec X,B,V,P,R,RP,T,T1,Q,U,D,AUQ;
23: X = ksp->vec_sol;
24: B = ksp->vec_rhs;
25: R = ksp->work[0];
26: RP = ksp->work[1];
27: V = ksp->work[2];
28: T = ksp->work[3];
29: Q = ksp->work[4];
30: P = ksp->work[5];
31: U = ksp->work[6];
32: D = ksp->work[7];
33: T1 = ksp->work[8];
34: AUQ = V;
36: /* Compute initial preconditioned residual */
37: KSPInitialResidual(ksp,X,V,T,R,B);
39: /* Test for nothing to do */
40: VecNorm(R,NORM_2,&dp);
41: PetscObjectSAWsTakeAccess((PetscObject)ksp);
42: ksp->rnorm = dp;
43: ksp->its = 0;
44: PetscObjectSAWsGrantAccess((PetscObject)ksp);
45: KSPMonitor(ksp,0,dp);
46: (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);
47: if (ksp->reason) return(0);
49: /* Make the initial Rp == R */
50: VecCopy(R,RP);
52: /* Set the initial conditions */
53: etaold = 0.0;
54: psiold = 0.0;
55: tau = dp;
56: dpold = dp;
58: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
59: VecCopy(R,U);
60: VecCopy(R,P);
61: KSP_PCApplyBAorAB(ksp,P,V,T);
62: VecSet(D,0.0);
64: i=0;
65: do {
66: PetscObjectSAWsTakeAccess((PetscObject)ksp);
67: ksp->its++;
68: PetscObjectSAWsGrantAccess((PetscObject)ksp);
69: VecDot(V,RP,&s); /* s <- (v,rp) */
70: a = rhoold / s; /* a <- rho / s */
71: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
72: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
73: KSP_PCApplyBAorAB(ksp,T,AUQ,T1);
74: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
75: VecNorm(R,NORM_2,&dp);
76: for (m=0; m<2; m++) {
77: if (!m) w = PetscSqrtReal(dp*dpold);
78: else w = dp;
79: psi = w / tau;
80: cm = 1.0 / PetscSqrtReal(1.0 + psi * psi);
81: tau = tau * psi * cm;
82: eta = cm * cm * a;
83: cf = psiold * psiold * etaold / a;
84: if (!m) {
85: VecAYPX(D,cf,U);
86: } else {
87: VecAYPX(D,cf,Q);
88: }
89: VecAXPY(X,eta,D);
91: dpest = PetscSqrtReal(m + 1.0) * tau;
92: PetscObjectSAWsTakeAccess((PetscObject)ksp);
93: ksp->rnorm = dpest;
94: PetscObjectSAWsGrantAccess((PetscObject)ksp);
95: KSPLogResidualHistory(ksp,dpest);
96: KSPMonitor(ksp,i+1,dpest);
97: (*ksp->converged)(ksp,i+1,dpest,&ksp->reason,ksp->cnvP);
98: if (ksp->reason) break;
100: etaold = eta;
101: psiold = psi;
102: }
103: if (ksp->reason) break;
105: VecDot(R,RP,&rho); /* rho <- (r,rp) */
106: b = rho / rhoold; /* b <- rho / rhoold */
107: VecWAXPY(U,b,Q,R); /* u <- r + b q */
108: VecAXPY(Q,b,P);
109: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
110: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
112: rhoold = rho;
113: dpold = dp;
115: i++;
116: } while (i<ksp->max_it);
117: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
119: KSPUnwindPreconditioner(ksp,X,T);
120: return(0);
121: }
123: /*MC
124: KSPTFQMR - A transpose free QMR (quasi minimal residual),
126: Options Database Keys:
127: . see KSPSolve()
129: Level: beginner
131: Notes: Supports left and right preconditioning, but not symmetric
133: The "residual norm" computed in this algorithm is actually just an upper bound on the actual residual norm.
134: That is for left preconditioning it is a bound on the preconditioned residual and for right preconditioning
135: it is a bound on the true residual.
137: References:
138: . 1. - Freund, 1993
140: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTCQMR
141: M*/
142: PETSC_EXTERN PetscErrorCode KSPCreate_TFQMR(KSP ksp)
143: {
147: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
148: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
150: ksp->data = (void*)0;
151: ksp->ops->setup = KSPSetUp_TFQMR;
152: ksp->ops->solve = KSPSolve_TFQMR;
153: ksp->ops->destroy = KSPDestroyDefault;
154: ksp->ops->buildsolution = KSPBuildSolutionDefault;
155: ksp->ops->buildresidual = KSPBuildResidualDefault;
156: ksp->ops->setfromoptions = 0;
157: ksp->ops->view = 0;
158: return(0);
159: }