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 static PetscErrorCode KSPSetUp_FBCGSR(KSP ksp) 13 { 14 PetscErrorCode ierr; 15 16 PetscFunctionBegin; 17 ierr = KSPSetWorkVecs(ksp,8);CHKERRQ(ierr); 18 PetscFunctionReturn(0); 19 } 20 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 36 PetscFunctionBegin; 37 PetscCheckFalse(!ksp->vec_rhs->petscnative,PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Only coded for PETSc vectors"); 38 ierr = VecGetLocalSize(ksp->vec_sol,&N);CHKERRQ(ierr); 39 40 X = ksp->vec_sol; 41 B = ksp->vec_rhs; 42 P2 = ksp->work[0]; 43 44 /* The followings are involved in modified inner product calculations and vector updates */ 45 RP = ksp->work[1]; ierr = VecGetArray(RP,(PetscScalar**)&rp);CHKERRQ(ierr); ierr = VecRestoreArray(RP,NULL);CHKERRQ(ierr); 46 R = ksp->work[2]; ierr = VecGetArray(R,(PetscScalar**)&r);CHKERRQ(ierr); ierr = VecRestoreArray(R,NULL);CHKERRQ(ierr); 47 P = ksp->work[3]; ierr = VecGetArray(P,(PetscScalar**)&p);CHKERRQ(ierr); ierr = VecRestoreArray(P,NULL);CHKERRQ(ierr); 48 V = ksp->work[4]; ierr = VecGetArray(V,(PetscScalar**)&v);CHKERRQ(ierr); ierr = VecRestoreArray(V,NULL);CHKERRQ(ierr); 49 S = ksp->work[5]; ierr = VecGetArray(S,(PetscScalar**)&s);CHKERRQ(ierr); ierr = VecRestoreArray(S,NULL);CHKERRQ(ierr); 50 T = ksp->work[6]; ierr = VecGetArray(T,(PetscScalar**)&t);CHKERRQ(ierr); ierr = VecRestoreArray(T,NULL);CHKERRQ(ierr); 51 S2 = ksp->work[7]; ierr = VecGetArray(S2,(PetscScalar**)&s2);CHKERRQ(ierr); ierr = VecRestoreArray(S2,NULL);CHKERRQ(ierr); 52 53 /* Only supports right preconditioning */ 54 PetscCheckFalse(ksp->pc_side != PC_RIGHT,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 ierr = VecDuplicate(X,&bcgs->guess);CHKERRQ(ierr); 58 } 59 ierr = VecCopy(X,bcgs->guess);CHKERRQ(ierr); 60 } else { 61 ierr = VecSet(X,0.0);CHKERRQ(ierr); 62 } 63 64 /* Compute initial residual */ 65 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 66 ierr = PCSetUp(pc);CHKERRQ(ierr); 67 ierr = PCGetOperators(pc,&mat,NULL);CHKERRQ(ierr); 68 if (!ksp->guess_zero) { 69 ierr = KSP_MatMult(ksp,mat,X,P2);CHKERRQ(ierr); /* P2 is used as temporary storage */ 70 ierr = VecCopy(B,R);CHKERRQ(ierr); 71 ierr = VecAXPY(R,-1.0,P2);CHKERRQ(ierr); 72 } else { 73 ierr = VecCopy(B,R);CHKERRQ(ierr); 74 } 75 76 /* Test for nothing to do */ 77 ierr = VecNorm(R,NORM_2,&rho);CHKERRQ(ierr); 78 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 79 ksp->its = 0; 80 if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = rho; 81 else ksp->rnorm = 0; 82 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 83 ierr = KSPLogResidualHistory(ksp,ksp->rnorm);CHKERRQ(ierr); 84 ierr = KSPMonitor(ksp,0,ksp->rnorm);CHKERRQ(ierr); 85 ierr = (*ksp->converged)(ksp,0,ksp->rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 86 if (ksp->reason) PetscFunctionReturn(0); 87 88 /* Initialize iterates */ 89 ierr = VecCopy(R,RP);CHKERRQ(ierr); /* rp <- r */ 90 ierr = VecCopy(R,P);CHKERRQ(ierr); /* p <- r */ 91 92 /* Big loop */ 93 for (i=0; i<ksp->max_it; i++) { 94 95 /* matmult and pc */ 96 ierr = KSP_PCApply(ksp,P,P2);CHKERRQ(ierr); /* p2 <- K p */ 97 ierr = KSP_MatMult(ksp,mat,P2,V);CHKERRQ(ierr); /* v <- A p2 */ 98 99 /* inner prodcuts */ 100 if (i==0) { 101 tau = rho*rho; 102 ierr = VecDot(V,RP,&sigma);CHKERRQ(ierr); /* sigma <- (v,rp) */ 103 } else { 104 ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 105 tau = sigma = 0.0; 106 for (j=0; j<N; j++) { 107 tau += r[j]*rp[j]; /* tau <- (r,rp) */ 108 sigma += v[j]*rp[j]; /* sigma <- (v,rp) */ 109 } 110 ierr = PetscLogFlops(4.0*N);CHKERRQ(ierr); 111 ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 112 insums[0] = tau; 113 insums[1] = sigma; 114 ierr = PetscLogEventBegin(VEC_ReduceCommunication,0,0,0,0);CHKERRQ(ierr); 115 ierr = MPIU_Allreduce(insums,outsums,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));CHKERRMPI(ierr); 116 ierr = PetscLogEventEnd(VEC_ReduceCommunication,0,0,0,0);CHKERRQ(ierr); 117 tau = outsums[0]; 118 sigma = outsums[1]; 119 } 120 121 /* scalar update */ 122 alpha = tau / sigma; 123 124 /* vector update */ 125 ierr = VecWAXPY(S,-alpha,V,R);CHKERRQ(ierr); /* s <- r - alpha v */ 126 127 /* matmult and pc */ 128 ierr = KSP_PCApply(ksp,S,S2);CHKERRQ(ierr); /* s2 <- K s */ 129 ierr = KSP_MatMult(ksp,mat,S2,T);CHKERRQ(ierr); /* t <- A s2 */ 130 131 /* inner prodcuts */ 132 ierr = PetscLogEventBegin(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 133 xi1 = xi2 = xi3 = xi4 = 0.0; 134 for (j=0; j<N; j++) { 135 xi1 += s[j]*s[j]; /* xi1 <- (s,s) */ 136 xi2 += t[j]*s[j]; /* xi2 <- (t,s) */ 137 xi3 += t[j]*t[j]; /* xi3 <- (t,t) */ 138 xi4 += t[j]*rp[j]; /* xi4 <- (t,rp) */ 139 } 140 ierr = PetscLogFlops(8.0*N);CHKERRQ(ierr); 141 ierr = PetscLogEventEnd(VEC_ReduceArithmetic,0,0,0,0);CHKERRQ(ierr); 142 143 insums[0] = xi1; 144 insums[1] = xi2; 145 insums[2] = xi3; 146 insums[3] = xi4; 147 148 ierr = PetscLogEventBegin(VEC_ReduceCommunication,0,0,0,0);CHKERRQ(ierr); 149 ierr = MPIU_Allreduce(insums,outsums,4,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ksp));CHKERRMPI(ierr); 150 ierr = PetscLogEventEnd(VEC_ReduceCommunication,0,0,0,0);CHKERRQ(ierr); 151 xi1 = outsums[0]; 152 xi2 = outsums[1]; 153 xi3 = outsums[2]; 154 xi4 = outsums[3]; 155 156 /* test denominator */ 157 if ((xi3 == 0.0) || (sigma == 0.0)) { 158 PetscCheckFalse(ksp->errorifnotconverged,PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has failed due to zero inner product"); 159 else ksp->reason = KSP_DIVERGED_BREAKDOWN; 160 ierr = PetscInfo(ksp,"KSPSolve has failed due to zero inner product\n");CHKERRQ(ierr); 161 break; 162 } 163 164 /* scalar updates */ 165 omega = xi2 / xi3; 166 beta = -xi4 / sigma; 167 rho = PetscSqrtReal(PetscAbsScalar(xi1 - omega * xi2)); /* residual norm */ 168 169 /* vector updates */ 170 ierr = VecAXPBYPCZ(X,alpha,omega,1.0,P2,S2);CHKERRQ(ierr); /* x <- alpha * p2 + omega * s2 + x */ 171 172 /* convergence test */ 173 ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); 174 ksp->its++; 175 if (ksp->normtype != KSP_NORM_NONE) ksp->rnorm = rho; 176 else ksp->rnorm = 0; 177 ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); 178 ierr = KSPLogResidualHistory(ksp,ksp->rnorm);CHKERRQ(ierr); 179 ierr = KSPMonitor(ksp,i+1,ksp->rnorm);CHKERRQ(ierr); 180 ierr = (*ksp->converged)(ksp,i+1,ksp->rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); 181 if (ksp->reason) break; 182 183 /* vector updates */ 184 ierr = PetscLogEventBegin(VEC_Ops,0,0,0,0);CHKERRQ(ierr); 185 for (j=0; j<N; j++) { 186 r[j] = s[j] - omega * t[j]; /* r <- s - omega t */ 187 p[j] = r[j] + beta * (p[j] - omega * v[j]); /* p <- r + beta * (p - omega v) */ 188 } 189 ierr = PetscLogFlops(6.0*N);CHKERRQ(ierr); 190 ierr = PetscLogEventEnd(VEC_Ops,0,0,0,0);CHKERRQ(ierr); 191 192 } 193 194 if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; 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: 207 Only allow right preconditioning 208 209 .seealso: KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP, KSPBICG, KSPFBCGSL, KSPSetPCSide() 210 M*/ 211 PETSC_EXTERN PetscErrorCode KSPCreate_FBCGSR(KSP ksp) 212 { 213 PetscErrorCode ierr; 214 KSP_BCGS *bcgs; 215 216 PetscFunctionBegin; 217 ierr = PetscNewLog(ksp,&bcgs);CHKERRQ(ierr); 218 219 ksp->data = bcgs; 220 ksp->ops->setup = KSPSetUp_FBCGSR; 221 ksp->ops->solve = KSPSolve_FBCGSR; 222 ksp->ops->destroy = KSPDestroy_BCGS; 223 ksp->ops->reset = KSPReset_BCGS; 224 ksp->ops->buildsolution = KSPBuildSolution_BCGS; 225 ksp->ops->buildresidual = KSPBuildResidualDefault; 226 ksp->ops->setfromoptions = KSPSetFromOptions_BCGS; 227 ksp->pc_side = PC_RIGHT; /* set default PC side */ 228 229 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,3);CHKERRQ(ierr); 230 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);CHKERRQ(ierr); 231 ierr = KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,1);CHKERRQ(ierr); 232 PetscFunctionReturn(0); 233 } 234