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) SETERRQ1(((PetscObject)ksp)->comm,PETSC_ERR_SUP,"KSP fbcgsr does not support %s",PCSides[ksp->pc_side]); 58 if (!ksp->guess_zero) { 59 if (!bcgs->guess) { 60 ierr = VecDuplicate(X,&bcgs->guess);CHKERRQ(ierr); 61 } 62 ierr = VecCopy(X,bcgs->guess);CHKERRQ(ierr); 63 } else { 64 ierr = VecSet(X,0.0);CHKERRQ(ierr); 65 } 66 67 /* Compute initial residual */ 68 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 69 ierr = PCSetUp(pc);CHKERRQ(ierr); 70 if (!ksp->guess_zero) { 71 ierr = MatMult(pc->mat,X,P2);CHKERRQ(ierr); /* P2 is used as temporary storage */ 72 ierr = VecCopy(B,R);CHKERRQ(ierr); 73 ierr = VecAXPY(R,-1.0,P2);CHKERRQ(ierr); 74 } else { 75 ierr = VecCopy(B,R);CHKERRQ(ierr); 76 } 77 78 /* Test for nothing to do */ 79 if (ksp->normtype != KSP_NORM_NONE) { 80 ierr = VecNorm(R,NORM_2,&rho);CHKERRQ(ierr); 81 } 82 ierr = PetscObjectTakeAccess(ksp);CHKERRQ(ierr); 83 ksp->its = 0; 84 ksp->rnorm = rho; 85 ierr = PetscObjectGrantAccess(ksp);CHKERRQ(ierr); 86 KSPLogResidualHistory(ksp,rho); 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 = PCApply(pc,P,P2);CHKERRQ(ierr); /* p2 <- K p */ 100 ierr = MatMult(pc->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 PetscLogFlops(4.0*N); 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,((PetscObject)ksp)->comm);CHKERRQ(ierr); 118 ierr = MPI_Allreduce(insums,outsums,2,MPIU_SCALAR,MPIU_SUM,((PetscObject)ksp)->comm);CHKERRQ(ierr); 119 ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);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 = PCApply(pc,S,S2);CHKERRQ(ierr); /* s2 <- K s */ 132 ierr = MatMult(pc->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 PetscLogFlops(8.0*N); 144 ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 145 insums[0] = xi1; 146 insums[1] = xi2; 147 insums[2] = xi3; 148 insums[3] = xi4; 149 ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);CHKERRQ(ierr); 150 ierr = MPI_Allreduce(insums,outsums,4,MPIU_SCALAR,MPIU_SUM,((PetscObject)ksp)->comm);CHKERRQ(ierr); 151 ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,((PetscObject)ksp)->comm);CHKERRQ(ierr); 152 xi1 = outsums[0]; 153 xi2 = outsums[1]; 154 xi3 = outsums[2]; 155 xi4 = outsums[3]; 156 157 /* test denominator */ 158 if (xi3 == 0.0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Divide by zero"); 159 if (sigma == 0.0) SETERRQ(((PetscObject)ksp)->comm,PETSC_ERR_PLIB,"Divide by zero"); 160 161 /* scalar updates */ 162 omega = xi2 / xi3; 163 beta = - xi4 / sigma; 164 rho = PetscSqrtReal(PetscAbsScalar(xi1 - omega * xi2)); /* residual norm */ 165 166 /* vector updates */ 167 ierr = VecAXPBYPCZ(X,alpha,omega,1.0,P2,S2);CHKERRQ(ierr); /* x <- alpha * p2 + omega * s2 + x */ 168 169 /* convergence test */ 170 ierr = PetscObjectTakeAccess(ksp);CHKERRQ(ierr); 171 ksp->its++; 172 ksp->rnorm = rho; 173 ierr = PetscObjectGrantAccess(ksp);CHKERRQ(ierr); 174 KSPLogResidualHistory(ksp,rho); 175 ierr = KSPMonitor(ksp,i+1,rho);CHKERRQ(ierr); 176 ierr = (*ksp->converged)(ksp,i+1,rho,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 177 if (ksp->reason) break; 178 179 /* vector updates */ 180 ierr = PetscLogEventBegin(VEC_Ops,0,0,0,0);CHKERRQ(ierr); 181 for (j=0; j<N; j++) { 182 r[j] = s[j] - omega * t[j]; /* r <- s - omega t */ 183 p[j] = r[j] + beta * (p[j] - omega * v[j]); /* p <- r + beta * (p - omega v) */ 184 } 185 PetscLogFlops(6.0*N); 186 ierr = PetscLogEventEnd(VEC_Ops,0,0,0,0);CHKERRQ(ierr); 187 188 } 189 190 if (i >= ksp->max_it) { 191 ksp->reason = KSP_DIVERGED_ITS; 192 } 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 EXTERN_C_BEGIN 209 #undef __FUNCT__ 210 #define __FUNCT__ "KSPCreate_FBCGSR" 211 PetscErrorCode KSPCreate_FBCGSR(KSP ksp) 212 { 213 PetscErrorCode ierr; 214 KSP_BCGS *bcgs; 215 216 PetscFunctionBegin; 217 ierr = PetscNewLog(ksp,KSP_BCGS,&bcgs);CHKERRQ(ierr); 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 = KSPDefaultBuildResidual; 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,2);CHKERRQ(ierr); 229 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,1);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232 EXTERN_C_END 233