Actual source code: tfqmr.c
petsc-3.14.3 2021-01-09
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: KSPCheckNorm(ksp,dp);
42: PetscObjectSAWsTakeAccess((PetscObject)ksp);
43: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = dp;
44: else ksp->rnorm = 0.0;
45: ksp->its = 0;
46: PetscObjectSAWsGrantAccess((PetscObject)ksp);
47: KSPMonitor(ksp,0,ksp->rnorm);
48: (*ksp->converged)(ksp,0,ksp->rnorm,&ksp->reason,ksp->cnvP);
49: if (ksp->reason) return(0);
51: /* Make the initial Rp == R */
52: VecCopy(R,RP);
54: /* Set the initial conditions */
55: etaold = 0.0;
56: psiold = 0.0;
57: tau = dp;
58: dpold = dp;
60: VecDot(R,RP,&rhoold); /* rhoold = (r,rp) */
61: VecCopy(R,U);
62: VecCopy(R,P);
63: KSP_PCApplyBAorAB(ksp,P,V,T);
64: VecSet(D,0.0);
66: i=0;
67: do {
68: PetscObjectSAWsTakeAccess((PetscObject)ksp);
69: ksp->its++;
70: PetscObjectSAWsGrantAccess((PetscObject)ksp);
71: VecDot(V,RP,&s); /* s <- (v,rp) */
72: KSPCheckDot(ksp,s);
73: a = rhoold / s; /* a <- rho / s */
74: VecWAXPY(Q,-a,V,U); /* q <- u - a v */
75: VecWAXPY(T,1.0,U,Q); /* t <- u + q */
76: KSP_PCApplyBAorAB(ksp,T,AUQ,T1);
77: VecAXPY(R,-a,AUQ); /* r <- r - a K (u + q) */
78: VecNorm(R,NORM_2,&dp);
79: KSPCheckNorm(ksp,dp);
80: for (m=0; m<2; m++) {
81: if (!m) w = PetscSqrtReal(dp*dpold);
82: else w = dp;
83: psi = w / tau;
84: cm = 1.0 / PetscSqrtReal(1.0 + psi * psi);
85: tau = tau * psi * cm;
86: eta = cm * cm * a;
87: cf = psiold * psiold * etaold / a;
88: if (!m) {
89: VecAYPX(D,cf,U);
90: } else {
91: VecAYPX(D,cf,Q);
92: }
93: VecAXPY(X,eta,D);
95: dpest = PetscSqrtReal(m + 1.0) * tau;
96: PetscObjectSAWsTakeAccess((PetscObject)ksp);
97: if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = dpest;
98: else ksp->rnorm = 0.0;
99: PetscObjectSAWsGrantAccess((PetscObject)ksp);
100: KSPLogResidualHistory(ksp,ksp->rnorm);
101: KSPMonitor(ksp,i+1,ksp->rnorm);
102: (*ksp->converged)(ksp,i+1,ksp->rnorm,&ksp->reason,ksp->cnvP);
103: if (ksp->reason) break;
105: etaold = eta;
106: psiold = psi;
107: }
108: if (ksp->reason) break;
110: VecDot(R,RP,&rho); /* rho <- (r,rp) */
111: b = rho / rhoold; /* b <- rho / rhoold */
112: VecWAXPY(U,b,Q,R); /* u <- r + b q */
113: VecAXPY(Q,b,P);
114: VecWAXPY(P,b,Q,U); /* p <- u + b(q + b p) */
115: KSP_PCApplyBAorAB(ksp,P,V,Q); /* v <- K p */
117: rhoold = rho;
118: dpold = dp;
120: i++;
121: } while (i<ksp->max_it);
122: if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
124: KSPUnwindPreconditioner(ksp,X,T);
125: return(0);
126: }
128: /*MC
129: KSPTFQMR - A transpose free QMR (quasi minimal residual),
131: Options Database Keys:
132: . see KSPSolve()
134: Level: beginner
136: Notes:
137: Supports left and right preconditioning, but not symmetric
139: The "residual norm" computed in this algorithm is actually just an upper bound on the actual residual norm.
140: That is for left preconditioning it is a bound on the preconditioned residual and for right preconditioning
141: it is a bound on the true residual.
143: References:
144: . 1. - Freund, 1993
146: .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPTCQMR
147: M*/
148: PETSC_EXTERN PetscErrorCode KSPCreate_TFQMR(KSP ksp)
149: {
153: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
154: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
155: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,1);
156: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);
158: ksp->data = (void*)0;
159: ksp->ops->setup = KSPSetUp_TFQMR;
160: ksp->ops->solve = KSPSolve_TFQMR;
161: ksp->ops->destroy = KSPDestroyDefault;
162: ksp->ops->buildsolution = KSPBuildSolutionDefault;
163: ksp->ops->buildresidual = KSPBuildResidualDefault;
164: ksp->ops->setfromoptions = NULL;
165: ksp->ops->view = NULL;
166: return(0);
167: }