xref: /petsc/src/ksp/ksp/impls/bcgs/fbcgsr/fbcgsr.c (revision 9c334d8fdb557fc53fd345d68cbb3545b09ccab8)
1 
2 /*
3     This file implements FBiCGStab-R.
4     Only allow right preconditioning.
5     FBiCGStab-R is a mathematically equivalent variant of FBiCGStab. Differences are:
6       (1) There are fewer MPI_Allreduce calls.
7       (2) The convergence occasionally is much faster than that of FBiCGStab.
8 */
9 #include <../src/ksp/ksp/impls/bcgs/bcgsimpl.h>       /*I  "petscksp.h"  I*/
10 #include <petsc/private/vecimpl.h>
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "KSPSetUp_FBCGSR"
14 static PetscErrorCode KSPSetUp_FBCGSR(KSP ksp)
15 {
16   PetscErrorCode ierr;
17 
18   PetscFunctionBegin;
19   ierr = KSPSetWorkVecs(ksp,8);CHKERRQ(ierr);
20   PetscFunctionReturn(0);
21 }
22 
23 #include <petsc/private/pcimpl.h>            /*I "petscksp.h" I*/
24 #undef __FUNCT__
25 #define __FUNCT__ "KSPSolve_FBCGSR"
26 static PetscErrorCode  KSPSolve_FBCGSR(KSP ksp)
27 {
28   PetscErrorCode    ierr;
29   PetscInt          i,j,N;
30   PetscScalar       tau,sigma,alpha,omega,beta;
31   PetscReal         rho;
32   PetscScalar       xi1,xi2,xi3,xi4;
33   Vec               X,B,P,P2,RP,R,V,S,T,S2;
34   PetscScalar       *PETSC_RESTRICT rp, *PETSC_RESTRICT r, *PETSC_RESTRICT p;
35   PetscScalar       *PETSC_RESTRICT v, *PETSC_RESTRICT s, *PETSC_RESTRICT t, *PETSC_RESTRICT s2;
36   PetscScalar       insums[4],outsums[4];
37   KSP_BCGS          *bcgs = (KSP_BCGS*)ksp->data;
38   PC                pc;
39 
40   PetscFunctionBegin;
41   if (!ksp->vec_rhs->petscnative) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Only coded for PETSc vectors");
42   ierr = VecGetLocalSize(ksp->vec_sol,&N);CHKERRQ(ierr);
43 
44   X  = ksp->vec_sol;
45   B  = ksp->vec_rhs;
46   P2 = ksp->work[0];
47 
48   /* The followings are involved in modified inner product calculations and vector updates */
49   RP = ksp->work[1]; ierr = VecGetArray(RP,(PetscScalar**)&rp);CHKERRQ(ierr); ierr = VecRestoreArray(RP,NULL);CHKERRQ(ierr);
50   R  = ksp->work[2]; ierr = VecGetArray(R,(PetscScalar**)&r);CHKERRQ(ierr);   ierr = VecRestoreArray(R,NULL);CHKERRQ(ierr);
51   P  = ksp->work[3]; ierr = VecGetArray(P,(PetscScalar**)&p);CHKERRQ(ierr);   ierr = VecRestoreArray(P,NULL);CHKERRQ(ierr);
52   V  = ksp->work[4]; ierr = VecGetArray(V,(PetscScalar**)&v);CHKERRQ(ierr);   ierr = VecRestoreArray(V,NULL);CHKERRQ(ierr);
53   S  = ksp->work[5]; ierr = VecGetArray(S,(PetscScalar**)&s);CHKERRQ(ierr);   ierr = VecRestoreArray(S,NULL);CHKERRQ(ierr);
54   T  = ksp->work[6]; ierr = VecGetArray(T,(PetscScalar**)&t);CHKERRQ(ierr);   ierr = VecRestoreArray(T,NULL);CHKERRQ(ierr);
55   S2 = ksp->work[7]; ierr = VecGetArray(S2,(PetscScalar**)&s2);CHKERRQ(ierr); ierr = VecRestoreArray(S2,NULL);CHKERRQ(ierr);
56 
57   /* Only supports right preconditioning */
58   if (ksp->pc_side != PC_RIGHT) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSP fbcgsr does not support %s",PCSides[ksp->pc_side]);
59   if (!ksp->guess_zero) {
60     if (!bcgs->guess) {
61       ierr = VecDuplicate(X,&bcgs->guess);CHKERRQ(ierr);
62     }
63     ierr = VecCopy(X,bcgs->guess);CHKERRQ(ierr);
64   } else {
65     ierr = VecSet(X,0.0);CHKERRQ(ierr);
66   }
67 
68   /* Compute initial residual */
69   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
70   ierr = PCSetUp(pc);CHKERRQ(ierr);
71   if (!ksp->guess_zero) {
72     ierr = KSP_MatMult(ksp,pc->mat,X,P2);CHKERRQ(ierr); /* P2 is used as temporary storage */
73     ierr = VecCopy(B,R);CHKERRQ(ierr);
74     ierr = VecAXPY(R,-1.0,P2);CHKERRQ(ierr);
75   } else {
76     ierr = VecCopy(B,R);CHKERRQ(ierr);
77   }
78 
79   /* Test for nothing to do */
80   ierr = VecNorm(R,NORM_2,&rho);CHKERRQ(ierr);
81   ierr       = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
82   ksp->its   = 0;
83   ksp->rnorm = rho;
84   ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
85   ierr = KSPLogResidualHistory(ksp,rho);CHKERRQ(ierr);
86   ierr = KSPMonitor(ksp,0,rho);CHKERRQ(ierr);
87   ierr = (*ksp->converged)(ksp,0,rho,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
88   if (ksp->reason) PetscFunctionReturn(0);
89 
90   /* Initialize iterates */
91   ierr = VecCopy(R,RP);CHKERRQ(ierr); /* rp <- r */
92   ierr = VecCopy(R,P);CHKERRQ(ierr); /* p <- r */
93 
94   /* Big loop */
95   for (i=0; i<ksp->max_it; i++) {
96 
97     /* matmult and pc */
98     ierr = KSP_PCApply(ksp,P,P2);CHKERRQ(ierr); /* p2 <- K p */
99     ierr = KSP_MatMult(ksp,pc->mat,P2,V);CHKERRQ(ierr); /* v <- A p2 */
100 
101     /* inner prodcuts */
102     if (i==0) {
103       tau  = rho*rho;
104       ierr = VecDot(V,RP,&sigma);CHKERRQ(ierr); /* sigma <- (v,rp) */
105     } else {
106       ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
107       tau  = sigma = 0.0;
108       for (j=0; j<N; j++) {
109         tau   += r[j]*rp[j]; /* tau <- (r,rp) */
110         sigma += v[j]*rp[j]; /* sigma <- (v,rp) */
111       }
112       ierr = PetscLogFlops(4.0*N);CHKERRQ(ierr);
113       ierr      = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
114       insums[0] = tau;
115       insums[1] = sigma;
116       ierr      = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr);
117       ierr      = MPIU_Allreduce(insums,outsums,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr);
118       ierr      = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr);
119       tau       = outsums[0];
120       sigma     = outsums[1];
121     }
122 
123     /* scalar update */
124     alpha = tau / sigma;
125 
126     /* vector update */
127     ierr = VecWAXPY(S,-alpha,V,R);CHKERRQ(ierr);  /* s <- r - alpha v */
128 
129     /* matmult and pc */
130     ierr = KSP_PCApply(ksp,S,S2);CHKERRQ(ierr); /* s2 <- K s */
131     ierr = KSP_MatMult(ksp,pc->mat,S2,T);CHKERRQ(ierr); /* t <- A s2 */
132 
133     /* inner prodcuts */
134     ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
135     xi1  = xi2 = xi3 = xi4 = 0.0;
136     for (j=0; j<N; j++) {
137       xi1 += s[j]*s[j]; /* xi1 <- (s,s) */
138       xi2 += t[j]*s[j]; /* xi2 <- (t,s) */
139       xi3 += t[j]*t[j]; /* xi3 <- (t,t) */
140       xi4 += t[j]*rp[j]; /* xi4 <- (t,rp) */
141     }
142     ierr = PetscLogFlops(8.0*N);CHKERRQ(ierr);
143     ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
144 
145     insums[0] = xi1;
146     insums[1] = xi2;
147     insums[2] = xi3;
148     insums[3] = xi4;
149 
150     ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr);
151     ierr = MPIU_Allreduce(insums,outsums,4,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr);
152     ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr);
153     xi1  = outsums[0];
154     xi2  = outsums[1];
155     xi3  = outsums[2];
156     xi4  = outsums[3];
157 
158     /* test denominator */
159     if (xi3 == 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Divide by zero");
160     if (sigma == 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Divide by zero");
161 
162     /* scalar updates */
163     omega = xi2 / xi3;
164     beta  = -xi4 / sigma;
165     rho   = PetscSqrtReal(PetscAbsScalar(xi1 - omega * xi2)); /* residual norm */
166 
167     /* vector updates */
168     ierr = VecAXPBYPCZ(X,alpha,omega,1.0,P2,S2);CHKERRQ(ierr); /* x <- alpha * p2 + omega * s2 + x */
169 
170     /* convergence test */
171     ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr);
172     ksp->its++;
173     ksp->rnorm = rho;
174     ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr);
175     ierr = KSPLogResidualHistory(ksp,rho);CHKERRQ(ierr);
176     ierr = KSPMonitor(ksp,i+1,rho);CHKERRQ(ierr);
177     ierr = (*ksp->converged)(ksp,i+1,rho,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
178     if (ksp->reason) break;
179 
180     /* vector updates */
181     ierr = PetscLogEventBegin(VEC_Ops,0,0,0,0);CHKERRQ(ierr);
182     for (j=0; j<N; j++) {
183       r[j] = s[j] - omega * t[j]; /* r <- s - omega t */
184       p[j] = r[j] + beta * (p[j] - omega * v[j]); /* p <- r + beta * (p - omega v) */
185     }
186     ierr = PetscLogFlops(6.0*N);CHKERRQ(ierr);
187     ierr = PetscLogEventEnd(VEC_Ops,0,0,0,0);CHKERRQ(ierr);
188 
189   }
190 
191   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
192   PetscFunctionReturn(0);
193 }
194 
195 /*MC
196      KSPFBCGSR - Implements a mathematically equivalent variant of FBiCGSTab.
197 
198    Options Database Keys:
199 .   see KSPSolve()
200 
201    Level: beginner
202 
203    Notes: Only allow right preconditioning
204 
205 .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPFBCGSL, KSPSetPCSide()
206 M*/
207 #undef __FUNCT__
208 #define __FUNCT__ "KSPCreate_FBCGSR"
209 PETSC_EXTERN PetscErrorCode KSPCreate_FBCGSR(KSP ksp)
210 {
211   PetscErrorCode ierr;
212   KSP_BCGS       *bcgs;
213 
214   PetscFunctionBegin;
215   ierr = PetscNewLog(ksp,&bcgs);CHKERRQ(ierr);
216 
217   ksp->data                = bcgs;
218   ksp->ops->setup          = KSPSetUp_FBCGSR;
219   ksp->ops->solve          = KSPSolve_FBCGSR;
220   ksp->ops->destroy        = KSPDestroy_BCGS;
221   ksp->ops->reset          = KSPReset_BCGS;
222   ksp->ops->buildsolution  = KSPBuildSolution_BCGS;
223   ksp->ops->buildresidual  = KSPBuildResidualDefault;
224   ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
225   ksp->pc_side             = PC_RIGHT; /* set default PC side */
226 
227   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr);
228   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);CHKERRQ(ierr);
229   PetscFunctionReturn(0);
230 }
231