xref: /petsc/src/ksp/ksp/impls/bcgs/fbcgsr/fbcgsr.c (revision deb2cd25b4974e90da8065ee93cb8b2449f9ae06)
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 
11 #undef __FUNCT__
12 #define __FUNCT__ "KSPSetUp_FBCGSR"
13 PetscErrorCode KSPSetUp_FBCGSR(KSP ksp)
14 {
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   ierr = KSPDefaultGetWork(ksp,8);CHKERRQ(ierr);
19   PetscFunctionReturn(0);
20 }
21 
22 #include <petsc-private/pcimpl.h>            /*I "petscksp.h" I*/
23 #undef __FUNCT__
24 #define __FUNCT__ "KSPSolve_FBCGSR"
25 PetscErrorCode  KSPSolve_FBCGSR(KSP ksp)
26 {
27   PetscErrorCode    ierr;
28   PetscInt          i,j,N;
29   PetscScalar       tau,sigma,alpha,omega,beta;
30   PetscReal         rho;
31   PetscScalar       xi1,xi2,xi3,xi4;
32   Vec               X,B,P,P2,RP,R,V,S,T,S2;
33   PetscScalar       *PETSC_RESTRICT rp, *PETSC_RESTRICT r, *PETSC_RESTRICT p;
34   PetscScalar       *PETSC_RESTRICT v, *PETSC_RESTRICT s, *PETSC_RESTRICT t, *PETSC_RESTRICT s2;
35   PetscScalar       insums[4],outsums[4];
36   KSP_BCGS          *bcgs = (KSP_BCGS*)ksp->data;
37   PC                pc;
38 
39   PetscFunctionBegin;
40   if (!ksp->vec_rhs->petscnative) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"Only coded for PETSc vectors");
41   ierr = VecGetLocalSize(ksp->vec_sol,&N);CHKERRQ(ierr);
42 
43   X   = ksp->vec_sol;
44   B   = ksp->vec_rhs;
45   P2  = ksp->work[0];
46 
47   /* The followings are involved in modified inner product calculations and vector updates */
48   RP  = ksp->work[1]; ierr = VecGetArray(RP,(PetscScalar**)&rp);CHKERRQ(ierr); ierr = VecRestoreArray(RP,PETSC_NULL);CHKERRQ(ierr);
49   R   = ksp->work[2]; ierr = VecGetArray(R,(PetscScalar**)&r);CHKERRQ(ierr);   ierr = VecRestoreArray(R,PETSC_NULL);CHKERRQ(ierr);
50   P   = ksp->work[3]; ierr = VecGetArray(P,(PetscScalar**)&p);CHKERRQ(ierr);   ierr = VecRestoreArray(P,PETSC_NULL);CHKERRQ(ierr);
51   V   = ksp->work[4]; ierr = VecGetArray(V,(PetscScalar**)&v);CHKERRQ(ierr);   ierr = VecRestoreArray(V,PETSC_NULL);CHKERRQ(ierr);
52   S   = ksp->work[5]; ierr = VecGetArray(S,(PetscScalar**)&s);CHKERRQ(ierr);   ierr = VecRestoreArray(S,PETSC_NULL);CHKERRQ(ierr);
53   T   = ksp->work[6]; ierr = VecGetArray(T,(PetscScalar**)&t);CHKERRQ(ierr);   ierr = VecRestoreArray(T,PETSC_NULL);CHKERRQ(ierr);
54   S2  = ksp->work[7]; ierr = VecGetArray(S2,(PetscScalar**)&s2);CHKERRQ(ierr); ierr = VecRestoreArray(S2,PETSC_NULL);CHKERRQ(ierr);
55 
56   /* Only supports right preconditioning */
57   if (ksp->pc_side != PC_RIGHT)
58     SETERRQ1(((PetscObject)ksp)->comm,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 = MatMult(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   if (ksp->normtype != KSP_NORM_NONE) {
81     ierr = VecNorm(R,NORM_2,&rho);CHKERRQ(ierr);
82   }
83   ierr = PetscObjectTakeAccess(ksp);CHKERRQ(ierr);
84   ksp->its   = 0;
85   ksp->rnorm = rho;
86   ierr = PetscObjectGrantAccess(ksp);CHKERRQ(ierr);
87   KSPLogResidualHistory(ksp,rho);
88   ierr = KSPMonitor(ksp,0,rho);CHKERRQ(ierr);
89   ierr = (*ksp->converged)(ksp,0,rho,&ksp->reason,ksp->cnvP);CHKERRQ(ierr);
90   if (ksp->reason) PetscFunctionReturn(0);
91 
92   /* Initialize iterates */
93   ierr = VecCopy(R,RP);CHKERRQ(ierr); /* rp <- r */
94   ierr = VecCopy(R,P);CHKERRQ(ierr); /* p <- r */
95 
96   /* Big loop */
97   for (i=0; i<ksp->max_it; i++) {
98 
99     /* matmult and pc */
100     ierr = PCApply(pc,P,P2);CHKERRQ(ierr); /* p2 <- K p */
101     ierr = MatMult(pc->mat,P2,V);CHKERRQ(ierr); /* v <- A p2 */
102 
103     /* inner prodcuts */
104     if (i==0) {
105       tau = rho*rho;
106       ierr = VecDot(V,RP,&sigma);CHKERRQ(ierr); /* sigma <- (v,rp) */
107     } else {
108       ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
109       tau = sigma = 0.0;
110       for (j=0; j<N; j++) {
111         tau += r[j]*rp[j]; /* tau <- (r,rp) */
112         sigma += v[j]*rp[j]; /* sigma <- (v,rp) */
113       }
114       PetscLogFlops(4.0*N);
115       ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
116       insums[0] = tau;
117       insums[1] = sigma;
118       ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);CHKERRQ(ierr);
119       ierr = MPI_Allreduce(insums,outsums,2,MPIU_SCALAR,MPIU_SUM,((PetscObject)ksp)->comm);CHKERRQ(ierr);
120       ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);CHKERRQ(ierr);
121       tau = outsums[0];
122       sigma = outsums[1];
123     }
124 
125     /* scalar update */
126     alpha = tau / sigma;
127 
128     /* vector update */
129     ierr = VecWAXPY(S,-alpha,V,R);CHKERRQ(ierr);  /* s <- r - alpha v */
130 
131     /* matmult and pc */
132     ierr = PCApply(pc,S,S2);CHKERRQ(ierr); /* s2 <- K s */
133     ierr = MatMult(pc->mat,S2,T);CHKERRQ(ierr); /* t <- A s2 */
134 
135     /* inner prodcuts */
136     ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
137     xi1 = xi2 = xi3 = xi4 = 0.0;
138     for (j=0; j<N; j++) {
139       xi1 += s[j]*s[j]; /* xi1 <- (s,s) */
140       xi2 += t[j]*s[j]; /* xi2 <- (t,s) */
141       xi3 += t[j]*t[j]; /* xi3 <- (t,t) */
142       xi4 += t[j]*rp[j]; /* xi4 <- (t,rp) */
143     }
144     PetscLogFlops(8.0*N);
145     ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr);
146     insums[0] = xi1;
147     insums[1] = xi2;
148     insums[2] = xi3;
149     insums[3] = xi4;
150     ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);CHKERRQ(ierr);
151     ierr = MPI_Allreduce(insums,outsums,4,MPIU_SCALAR,MPIU_SUM,((PetscObject)ksp)->comm);CHKERRQ(ierr);
152     ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);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(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Divide by zero");
160     if (sigma == 0.0) SETERRQ(((PetscObject)ksp)->comm,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 = PetscObjectTakeAccess(ksp);CHKERRQ(ierr);
172     ksp->its++;
173     ksp->rnorm = rho;
174     ierr = PetscObjectGrantAccess(ksp);CHKERRQ(ierr);
175     KSPLogResidualHistory(ksp,rho);
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     PetscLogFlops(6.0*N);
187     ierr = PetscLogEventEnd(VEC_Ops,0,0,0,0);CHKERRQ(ierr);
188 
189   }
190 
191   if (i >= ksp->max_it) {
192     ksp->reason = KSP_DIVERGED_ITS;
193   }
194 
195   PetscFunctionReturn(0);
196 }
197 
198 /*MC
199      KSPFBCGSR - Implements a mathematically equivalent variant of FBiCGSTab.
200 
201    Options Database Keys:
202 .   see KSPSolve()
203 
204    Level: beginner
205 
206    Notes: Only allow right preconditioning
207 
208 .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPFBCGSL, KSPSetPCSide()
209 M*/
210 EXTERN_C_BEGIN
211 #undef __FUNCT__
212 #define __FUNCT__ "KSPCreate_FBCGSR"
213 PetscErrorCode  KSPCreate_FBCGSR(KSP ksp)
214 {
215   PetscErrorCode ierr;
216   KSP_BCGS       *bcgs;
217 
218   PetscFunctionBegin;
219   ierr = PetscNewLog(ksp,KSP_BCGS,&bcgs);CHKERRQ(ierr);
220   ksp->data                 = bcgs;
221   ksp->ops->setup           = KSPSetUp_FBCGSR;
222   ksp->ops->solve           = KSPSolve_FBCGSR;
223   ksp->ops->destroy         = KSPDestroy_BCGS;
224   ksp->ops->reset           = KSPReset_BCGS;
225   ksp->ops->buildsolution   = KSPBuildSolution_BCGS;
226   ksp->ops->buildresidual   = KSPDefaultBuildResidual;
227   ksp->ops->setfromoptions  = KSPSetFromOptions_BCGS;
228   ksp->ops->view            = KSPView_BCGS;
229   ksp->pc_side              = PC_RIGHT; /* set default PC side */
230 
231   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);CHKERRQ(ierr);
232   ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);CHKERRQ(ierr);
233   PetscFunctionReturn(0);
234 }
235 EXTERN_C_END
236