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 #undef __FUNCT__ 24 #define __FUNCT__ "KSPSolve_FBCGSR" 25 static 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 Mat mat; 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 ierr = PCGetOperators(pc,&mat,NULL);CHKERRQ(ierr); 72 if (!ksp->guess_zero) { 73 ierr = KSP_MatMult(ksp,mat,X,P2);CHKERRQ(ierr); /* P2 is used as temporary storage */ 74 ierr = VecCopy(B,R);CHKERRQ(ierr); 75 ierr = VecAXPY(R,-1.0,P2);CHKERRQ(ierr); 76 } else { 77 ierr = VecCopy(B,R);CHKERRQ(ierr); 78 } 79 80 /* Test for nothing to do */ 81 ierr = VecNorm(R,NORM_2,&rho);CHKERRQ(ierr); 82 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 83 ksp->its = 0; 84 ksp->rnorm = rho; 85 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 86 ierr = KSPLogResidualHistory(ksp,rho);CHKERRQ(ierr); 87 ierr = KSPMonitor(ksp,0,rho);CHKERRQ(ierr); 88 ierr = (*ksp->converged)(ksp,0,rho,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 89 if (ksp->reason) PetscFunctionReturn(0); 90 91 /* Initialize iterates */ 92 ierr = VecCopy(R,RP);CHKERRQ(ierr); /* rp <- r */ 93 ierr = VecCopy(R,P);CHKERRQ(ierr); /* p <- r */ 94 95 /* Big loop */ 96 for (i=0; i<ksp->max_it; i++) { 97 98 /* matmult and pc */ 99 ierr = KSP_PCApply(ksp,P,P2);CHKERRQ(ierr); /* p2 <- K p */ 100 ierr = KSP_MatMult(ksp,mat,P2,V);CHKERRQ(ierr); /* v <- A p2 */ 101 102 /* inner prodcuts */ 103 if (i==0) { 104 tau = rho*rho; 105 ierr = VecDot(V,RP,&sigma);CHKERRQ(ierr); /* sigma <- (v,rp) */ 106 } else { 107 ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 108 tau = sigma = 0.0; 109 for (j=0; j<N; j++) { 110 tau += r[j]*rp[j]; /* tau <- (r,rp) */ 111 sigma += v[j]*rp[j]; /* sigma <- (v,rp) */ 112 } 113 ierr = PetscLogFlops(4.0*N);CHKERRQ(ierr); 114 ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 115 insums[0] = tau; 116 insums[1] = sigma; 117 ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr); 118 ierr = MPIU_Allreduce(insums,outsums,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr); 119 ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr); 120 tau = outsums[0]; 121 sigma = outsums[1]; 122 } 123 124 /* scalar update */ 125 alpha = tau / sigma; 126 127 /* vector update */ 128 ierr = VecWAXPY(S,-alpha,V,R);CHKERRQ(ierr); /* s <- r - alpha v */ 129 130 /* matmult and pc */ 131 ierr = KSP_PCApply(ksp,S,S2);CHKERRQ(ierr); /* s2 <- K s */ 132 ierr = KSP_MatMult(ksp,mat,S2,T);CHKERRQ(ierr); /* t <- A s2 */ 133 134 /* inner prodcuts */ 135 ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 136 xi1 = xi2 = xi3 = xi4 = 0.0; 137 for (j=0; j<N; j++) { 138 xi1 += s[j]*s[j]; /* xi1 <- (s,s) */ 139 xi2 += t[j]*s[j]; /* xi2 <- (t,s) */ 140 xi3 += t[j]*t[j]; /* xi3 <- (t,t) */ 141 xi4 += t[j]*rp[j]; /* xi4 <- (t,rp) */ 142 } 143 ierr = PetscLogFlops(8.0*N);CHKERRQ(ierr); 144 ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 145 146 insums[0] = xi1; 147 insums[1] = xi2; 148 insums[2] = xi3; 149 insums[3] = xi4; 150 151 ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr); 152 ierr = MPIU_Allreduce(insums,outsums,4,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr); 153 ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,PetscObjectComm((PetscObject)ksp));CHKERRQ(ierr); 154 xi1 = outsums[0]; 155 xi2 = outsums[1]; 156 xi3 = outsums[2]; 157 xi4 = outsums[3]; 158 159 /* test denominator */ 160 if (xi3 == 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Divide by zero"); 161 if (sigma == 0.0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Divide by zero"); 162 163 /* scalar updates */ 164 omega = xi2 / xi3; 165 beta = -xi4 / sigma; 166 rho = PetscSqrtReal(PetscAbsScalar(xi1 - omega * xi2)); /* residual norm */ 167 168 /* vector updates */ 169 ierr = VecAXPBYPCZ(X,alpha,omega,1.0,P2,S2);CHKERRQ(ierr); /* x <- alpha * p2 + omega * s2 + x */ 170 171 /* convergence test */ 172 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 173 ksp->its++; 174 ksp->rnorm = rho; 175 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 176 ierr = KSPLogResidualHistory(ksp,rho);CHKERRQ(ierr); 177 ierr = KSPMonitor(ksp,i+1,rho);CHKERRQ(ierr); 178 ierr = (*ksp->converged)(ksp,i+1,rho,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 179 if (ksp->reason) break; 180 181 /* vector updates */ 182 ierr = PetscLogEventBegin(VEC_Ops,0,0,0,0);CHKERRQ(ierr); 183 for (j=0; j<N; j++) { 184 r[j] = s[j] - omega * t[j]; /* r <- s - omega t */ 185 p[j] = r[j] + beta * (p[j] - omega * v[j]); /* p <- r + beta * (p - omega v) */ 186 } 187 ierr = PetscLogFlops(6.0*N);CHKERRQ(ierr); 188 ierr = PetscLogEventEnd(VEC_Ops,0,0,0,0);CHKERRQ(ierr); 189 190 } 191 192 if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; 193 PetscFunctionReturn(0); 194 } 195 196 /*MC 197 KSPFBCGSR - Implements a mathematically equivalent variant of FBiCGSTab. 198 199 Options Database Keys: 200 . see KSPSolve() 201 202 Level: beginner 203 204 Notes: Only allow right preconditioning 205 206 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPFBCGSL, KSPSetPCSide() 207 M*/ 208 #undef __FUNCT__ 209 #define __FUNCT__ "KSPCreate_FBCGSR" 210 PETSC_EXTERN PetscErrorCode KSPCreate_FBCGSR(KSP ksp) 211 { 212 PetscErrorCode ierr; 213 KSP_BCGS *bcgs; 214 215 PetscFunctionBegin; 216 ierr = PetscNewLog(ksp,&bcgs);CHKERRQ(ierr); 217 218 ksp->data = bcgs; 219 ksp->ops->setup = KSPSetUp_FBCGSR; 220 ksp->ops->solve = KSPSolve_FBCGSR; 221 ksp->ops->destroy = KSPDestroy_BCGS; 222 ksp->ops->reset = KSPReset_BCGS; 223 ksp->ops->buildsolution = KSPBuildSolution_BCGS; 224 ksp->ops->buildresidual = KSPBuildResidualDefault; 225 ksp->ops->setfromoptions = KSPSetFromOptions_BCGS; 226 ksp->pc_side = PC_RIGHT; /* set default PC side */ 227 228 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr); 229 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232