Actual source code: fbcgsr.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  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>
 10:  #include <petsc/private/vecimpl.h>

 12: static PetscErrorCode KSPSetUp_FBCGSR(KSP ksp)
 13: {

 17:   KSPSetWorkVecs(ksp,8);
 18:   return(0);
 19: }

 21: static PetscErrorCode  KSPSolve_FBCGSR(KSP ksp)
 22: {
 23:   PetscErrorCode    ierr;
 24:   PetscInt          i,j,N;
 25:   PetscScalar       tau,sigma,alpha,omega,beta;
 26:   PetscReal         rho;
 27:   PetscScalar       xi1,xi2,xi3,xi4;
 28:   Vec               X,B,P,P2,RP,R,V,S,T,S2;
 29:   PetscScalar       *PETSC_RESTRICT rp, *PETSC_RESTRICT r, *PETSC_RESTRICT p;
 30:   PetscScalar       *PETSC_RESTRICT v, *PETSC_RESTRICT s, *PETSC_RESTRICT t, *PETSC_RESTRICT s2;
 31:   PetscScalar       insums[4],outsums[4];
 32:   KSP_BCGS          *bcgs = (KSP_BCGS*)ksp->data;
 33:   PC                pc;
 34:   Mat               mat;
 35: 
 37:   if (!ksp->vec_rhs->petscnative) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Only coded for PETSc vectors");
 38:   VecGetLocalSize(ksp->vec_sol,&N);

 40:   X  = ksp->vec_sol;
 41:   B  = ksp->vec_rhs;
 42:   P2 = ksp->work[0];

 44:   /* The followings are involved in modified inner product calculations and vector updates */
 45:   RP = ksp->work[1]; VecGetArray(RP,(PetscScalar**)&rp); VecRestoreArray(RP,NULL);
 46:   R  = ksp->work[2]; VecGetArray(R,(PetscScalar**)&r);   VecRestoreArray(R,NULL);
 47:   P  = ksp->work[3]; VecGetArray(P,(PetscScalar**)&p);   VecRestoreArray(P,NULL);
 48:   V  = ksp->work[4]; VecGetArray(V,(PetscScalar**)&v);   VecRestoreArray(V,NULL);
 49:   S  = ksp->work[5]; VecGetArray(S,(PetscScalar**)&s);   VecRestoreArray(S,NULL);
 50:   T  = ksp->work[6]; VecGetArray(T,(PetscScalar**)&t);   VecRestoreArray(T,NULL);
 51:   S2 = ksp->work[7]; VecGetArray(S2,(PetscScalar**)&s2); VecRestoreArray(S2,NULL);

 53:   /* Only supports right preconditioning */
 54:   if (ksp->pc_side != PC_RIGHT) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"KSP fbcgsr does not support %s",PCSides[ksp->pc_side]);
 55:   if (!ksp->guess_zero) {
 56:     if (!bcgs->guess) {
 57:       VecDuplicate(X,&bcgs->guess);
 58:     }
 59:     VecCopy(X,bcgs->guess);
 60:   } else {
 61:     VecSet(X,0.0);
 62:   }

 64:   /* Compute initial residual */
 65:   KSPGetPC(ksp,&pc);
 66:   PCSetUp(pc);
 67:   PCGetOperators(pc,&mat,NULL);
 68:   if (!ksp->guess_zero) {
 69:     KSP_MatMult(ksp,mat,X,P2); /* P2 is used as temporary storage */
 70:     VecCopy(B,R);
 71:     VecAXPY(R,-1.0,P2);
 72:   } else {
 73:     VecCopy(B,R);
 74:   }

 76:   /* Test for nothing to do */
 77:   VecNorm(R,NORM_2,&rho);
 78:   PetscObjectSAWsTakeAccess((PetscObject)ksp);
 79:   ksp->its   = 0;
 80:   ksp->rnorm = rho;
 81:   PetscObjectSAWsGrantAccess((PetscObject)ksp);
 82:   KSPLogResidualHistory(ksp,rho);
 83:   KSPMonitor(ksp,0,rho);
 84:   (*ksp->converged)(ksp,0,rho,&ksp->reason,ksp->cnvP);
 85:   if (ksp->reason) return(0);

 87:   /* Initialize iterates */
 88:   VecCopy(R,RP); /* rp <- r */
 89:   VecCopy(R,P); /* p <- r */

 91:   /* Big loop */
 92:   for (i=0; i<ksp->max_it; i++) {

 94:     /* matmult and pc */
 95:     KSP_PCApply(ksp,P,P2); /* p2 <- K p */
 96:     KSP_MatMult(ksp,mat,P2,V); /* v <- A p2 */

 98:     /* inner prodcuts */
 99:     if (i==0) {
100:       tau  = rho*rho;
101:       VecDot(V,RP,&sigma); /* sigma <- (v,rp) */
102:     } else {
103:       PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);
104:       tau  = sigma = 0.0;
105:       for (j=0; j<N; j++) {
106:         tau   += r[j]*rp[j]; /* tau <- (r,rp) */
107:         sigma += v[j]*rp[j]; /* sigma <- (v,rp) */
108:       }
109:       PetscLogFlops(4.0*N);
110:       PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);
111:       insums[0] = tau;
112:       insums[1] = sigma;
113:       PetscLogEventBegin(VEC_ReduceCommunication,0,0,0,0);
114:       MPIU_Allreduce(insums,outsums,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));
115:       PetscLogEventEnd(VEC_ReduceCommunication,0,0,0,0);
116:       tau       = outsums[0];
117:       sigma     = outsums[1];
118:     }

120:     /* scalar update */
121:     alpha = tau / sigma;

123:     /* vector update */
124:     VecWAXPY(S,-alpha,V,R);  /* s <- r - alpha v */

126:     /* matmult and pc */
127:     KSP_PCApply(ksp,S,S2); /* s2 <- K s */
128:     KSP_MatMult(ksp,mat,S2,T); /* t <- A s2 */

130:     /* inner prodcuts */
131:     PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);
132:     xi1  = xi2 = xi3 = xi4 = 0.0;
133:     for (j=0; j<N; j++) {
134:       xi1 += s[j]*s[j]; /* xi1 <- (s,s) */
135:       xi2 += t[j]*s[j]; /* xi2 <- (t,s) */
136:       xi3 += t[j]*t[j]; /* xi3 <- (t,t) */
137:       xi4 += t[j]*rp[j]; /* xi4 <- (t,rp) */
138:     }
139:     PetscLogFlops(8.0*N);
140:     PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);

142:     insums[0] = xi1;
143:     insums[1] = xi2;
144:     insums[2] = xi3;
145:     insums[3] = xi4;

147:     PetscLogEventBegin(VEC_ReduceCommunication,0,0,0,0);
148:     MPIU_Allreduce(insums,outsums,4,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));
149:     PetscLogEventEnd(VEC_ReduceCommunication,0,0,0,0);
150:     xi1  = outsums[0];
151:     xi2  = outsums[1];
152:     xi3  = outsums[2];
153:     xi4  = outsums[3];

155:     /* test denominator */
156:     if (xi3 == 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Divide by zero");
157:     if (sigma == 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Divide by zero");

159:     /* scalar updates */
160:     omega = xi2 / xi3;
161:     beta  = -xi4 / sigma;
162:     rho   = PetscSqrtReal(PetscAbsScalar(xi1 - omega * xi2)); /* residual norm */

164:     /* vector updates */
165:     VecAXPBYPCZ(X,alpha,omega,1.0,P2,S2); /* x <- alpha * p2 + omega * s2 + x */

167:     /* convergence test */
168:     PetscObjectSAWsTakeAccess((PetscObject)ksp);
169:     ksp->its++;
170:     ksp->rnorm = rho;
171:     PetscObjectSAWsGrantAccess((PetscObject)ksp);
172:     KSPLogResidualHistory(ksp,rho);
173:     KSPMonitor(ksp,i+1,rho);
174:     (*ksp->converged)(ksp,i+1,rho,&ksp->reason,ksp->cnvP);
175:     if (ksp->reason) break;

177:     /* vector updates */
178:     PetscLogEventBegin(VEC_Ops,0,0,0,0);
179:     for (j=0; j<N; j++) {
180:       r[j] = s[j] - omega * t[j]; /* r <- s - omega t */
181:       p[j] = r[j] + beta * (p[j] - omega * v[j]); /* p <- r + beta * (p - omega v) */
182:     }
183:     PetscLogFlops(6.0*N);
184:     PetscLogEventEnd(VEC_Ops,0,0,0,0);

186:   }

188:   if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
189:   return(0);
190: }

192: /*MC
193:      KSPFBCGSR - Implements a mathematically equivalent variant of FBiCGSTab.

195:    Options Database Keys:
196: .   see KSPSolve()

198:    Level: beginner

200:    Notes:
201:     Only allow right preconditioning

203: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPFBCGSL, KSPSetPCSide()
204: M*/
205: PETSC_EXTERN PetscErrorCode KSPCreate_FBCGSR(KSP ksp)
206: {
208:   KSP_BCGS       *bcgs;

211:   PetscNewLog(ksp,&bcgs);

213:   ksp->data                = bcgs;
214:   ksp->ops->setup          = KSPSetUp_FBCGSR;
215:   ksp->ops->solve          = KSPSolve_FBCGSR;
216:   ksp->ops->destroy        = KSPDestroy_BCGS;
217:   ksp->ops->reset          = KSPReset_BCGS;
218:   ksp->ops->buildsolution  = KSPBuildSolution_BCGS;
219:   ksp->ops->buildresidual  = KSPBuildResidualDefault;
220:   ksp->ops->setfromoptions = KSPSetFromOptions_BCGS;
221:   ksp->pc_side             = PC_RIGHT; /* set default PC side */

223:   KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);
224:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
225:   return(0);
226: }