PETSC_EXTERN void kspconvergedskip_(KSP *ksp,PetscInt *n,PetscReal *rnorm,KSPConvergedReason *flag,void *dummy,PetscErrorCode *ierr) { CHKFORTRANNULLOBJECT(dummy); *ierr = KSPConvergedSkip(*ksp,*n,*rnorm,flag,dummy); }
static PetscErrorCode KSPSolve_PIPEFCG(KSP ksp) { PetscErrorCode ierr; KSP_PIPEFCG *pipefcg; PetscScalar gamma; PetscReal dp=0.0; Vec B,R,Z,X; Mat Amat,Pmat; #define VecXDot(x,y,a) (((pipefcg->type) == (KSP_CG_HERMITIAN)) ? VecDot (x,y,a) : VecTDot (x,y,a)) PetscFunctionBegin; ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); pipefcg = (KSP_PIPEFCG*)ksp->data; X = ksp->vec_sol; B = ksp->vec_rhs; R = ksp->work[0]; Z = ksp->work[1]; ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); /* Compute initial residual needed for convergence check*/ ksp->its = 0; if (!ksp->guess_zero) { ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr); ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr); /* r <- b - Ax */ } else { ierr = VecCopy(B,R);CHKERRQ(ierr); /* r <- b (x is 0) */ } switch (ksp->normtype) { case KSP_NORM_PRECONDITIONED: ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- dqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */ break; case KSP_NORM_UNPRECONDITIONED: ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e) */ break; case KSP_NORM_NATURAL: ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ ierr = VecXDot(Z,R,&gamma);CHKERRQ(ierr); dp = PetscSqrtReal(PetscAbsScalar(gamma)); /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e) */ break; case KSP_NORM_NONE: dp = 0.0; break; default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); } /* Initial Convergence Check */ ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); ierr = KSPMonitor(ksp,0,dp);CHKERRQ(ierr); ksp->rnorm = dp; if (ksp->normtype == KSP_NORM_NONE) { ierr = KSPConvergedSkip (ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); } else { ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); } if (ksp->reason) PetscFunctionReturn(0); do { /* A cycle is broken only if a norm breakdown occurs. If not the entire solve happens in a single cycle. This is coded this way to allow both truncation and truncation-restart strategy (see KSPFCDGetNumOldDirections()) */ ierr = KSPSolve_PIPEFCG_cycle(ksp);CHKERRQ(ierr); if (ksp->reason) break; if (pipefcg->norm_breakdown) { pipefcg->n_restarts++; pipefcg->norm_breakdown = PETSC_FALSE; } } while (ksp->its < ksp->max_it); if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; PetscFunctionReturn(0); }
PetscErrorCode KSPSolve_FCG(KSP ksp) { PetscErrorCode ierr; PetscInt i,k,idx,mi; KSP_FCG *fcg = (KSP_FCG*)ksp->data; PetscScalar alpha=0.0,beta = 0.0,dpi,s; PetscReal dp=0.0; Vec B,R,Z,X,Pcurr,Ccurr; Mat Amat,Pmat; PetscInt eigs = ksp->calc_sings; /* Variables for eigen estimation - START*/ PetscInt stored_max_it = ksp->max_it; PetscScalar alphaold = 0,betaold = 1.0,*e = 0,*d = 0;/* Variables for eigen estimation - FINISH */ PetscFunctionBegin; #define VecXDot(x,y,a) (((fcg->type) == (KSP_CG_HERMITIAN)) ? VecDot(x,y,a) : VecTDot(x,y,a)) #define VecXMDot(a,b,c,d) (((fcg->type) == (KSP_CG_HERMITIAN)) ? VecMDot(a,b,c,d) : VecMTDot(a,b,c,d)) X = ksp->vec_sol; B = ksp->vec_rhs; R = ksp->work[0]; Z = ksp->work[1]; ierr = PCGetOperators(ksp->pc,&Amat,&Pmat);CHKERRQ(ierr); if (eigs) {e = fcg->e; d = fcg->d; e[0] = 0.0; } /* Compute initial residual needed for convergence check*/ ksp->its = 0; if (!ksp->guess_zero) { ierr = KSP_MatMult(ksp,Amat,X,R);CHKERRQ(ierr); ierr = VecAYPX(R,-1.0,B);CHKERRQ(ierr); /* r <- b - Ax */ } else { ierr = VecCopy(B,R);CHKERRQ(ierr); /* r <- b (x is 0) */ } switch (ksp->normtype) { case KSP_NORM_PRECONDITIONED: ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- dqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */ break; case KSP_NORM_UNPRECONDITIONED: ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e) */ break; case KSP_NORM_NATURAL: ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ ierr = VecXDot(R,Z,&s);CHKERRQ(ierr); dp = PetscSqrtReal(PetscAbsScalar(s)); /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e) */ break; case KSP_NORM_NONE: dp = 0.0; break; default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); } /* Initial Convergence Check */ ierr = KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); ierr = KSPMonitor(ksp,0,dp);CHKERRQ(ierr); ksp->rnorm = dp; if (ksp->normtype == KSP_NORM_NONE) { ierr = KSPConvergedSkip(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); } else { ierr = (*ksp->converged)(ksp,0,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); } if (ksp->reason) PetscFunctionReturn(0); /* Apply PC if not already done for convergence check */ if(ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE){ ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ } i = 0; do { ksp->its = i+1; /* If needbe, allocate a new chunk of vectors in P and C */ ierr = KSPAllocateVectors_FCG(ksp,i+1,fcg->vecb);CHKERRQ(ierr); /* Note that we wrap around and start clobbering old vectors */ idx = i % (fcg->mmax+1); Pcurr = fcg->Pvecs[idx]; Ccurr = fcg->Cvecs[idx]; /* Compute a new column of P (Currently does not support modified G-S or iterative refinement)*/ switch(fcg->truncstrat){ case KSP_FCG_TRUNC_TYPE_NOTAY : mi = PetscMax(1,i%(fcg->mmax+1)); break; case KSP_FCG_TRUNC_TYPE_STANDARD : mi = fcg->mmax; break; default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unrecognized FCG Truncation Strategy");CHKERRQ(ierr); } ierr = VecCopy(Z,Pcurr);CHKERRQ(ierr); { PetscInt l,ndots; l = PetscMax(0,i-mi); ndots = i-l; if (ndots){ PetscInt j; Vec *Pold, *Cold; PetscScalar *dots; ierr = PetscMalloc3(ndots,&dots,ndots,&Cold,ndots,&Pold);CHKERRQ(ierr); for(k=l,j=0;j<ndots;++k,++j){ idx = k % (fcg->mmax+1); Cold[j] = fcg->Cvecs[idx]; Pold[j] = fcg->Pvecs[idx]; } ierr = VecXMDot(Z,ndots,Cold,dots);CHKERRQ(ierr); for(k=0;k<ndots;++k){ dots[k] = -dots[k]; } ierr = VecMAXPY(Pcurr,ndots,dots,Pold);CHKERRQ(ierr); ierr = PetscFree3(dots,Cold,Pold);CHKERRQ(ierr); } } /* Update X and R */ betaold = beta; ierr = VecXDot(Pcurr,R,&beta);CHKERRQ(ierr); /* beta <- pi'*r */ ierr = KSP_MatMult(ksp,Amat,Pcurr,Ccurr);CHKERRQ(ierr); /* w <- A*pi (stored in ci) */ ierr = VecXDot(Pcurr,Ccurr,&dpi);CHKERRQ(ierr); /* dpi <- pi'*w */ alphaold = alpha; alpha = beta / dpi; /* alpha <- beta/dpi */ ierr = VecAXPY(X,alpha,Pcurr);CHKERRQ(ierr); /* x <- x + alpha * pi */ ierr = VecAXPY(R,-alpha,Ccurr);CHKERRQ(ierr); /* r <- r - alpha * wi */ /* Compute norm for convergence check */ switch (ksp->normtype) { case KSP_NORM_PRECONDITIONED: ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ ierr = VecNorm(Z,NORM_2,&dp);CHKERRQ(ierr); /* dp <- sqrt(z'*z) = sqrt(e'*A'*B'*B*A*e) */ break; case KSP_NORM_UNPRECONDITIONED: ierr = VecNorm(R,NORM_2,&dp);CHKERRQ(ierr); /* dp <- sqrt(r'*r) = sqrt(e'*A'*A*e) */ break; case KSP_NORM_NATURAL: ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ ierr = VecXDot(R,Z,&s);CHKERRQ(ierr); dp = PetscSqrtReal(PetscAbsScalar(s)); /* dp <- sqrt(r'*z) = sqrt(e'*A'*B*A*e) */ break; case KSP_NORM_NONE: dp = 0.0; break; default: SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"%s",KSPNormTypes[ksp->normtype]); } /* Check for convergence */ ksp->rnorm = dp; KSPLogResidualHistory(ksp,dp);CHKERRQ(ierr); ierr = KSPMonitor(ksp,i+1,dp);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i+1,dp,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; /* Apply PC if not already done for convergence check */ if(ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->normtype == KSP_NORM_NONE){ ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /* z <- Br */ } /* Compute current C (which is W/dpi) */ ierr = VecScale(Ccurr,1.0/dpi);CHKERRQ(ierr); /* w <- ci/dpi */ if (eigs) { if (i > 0) { if (ksp->max_it != stored_max_it) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Can not change maxit AND calculate eigenvalues"); e[i] = PetscSqrtReal(PetscAbsScalar(beta/betaold))/alphaold; d[i] = PetscSqrtReal(PetscAbsScalar(beta/betaold))*e[i] + 1.0/alpha; } else { d[i] = PetscSqrtReal(PetscAbsScalar(beta))*e[i] + 1.0/alpha; } } ++i; } while (i<ksp->max_it); if (i >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS; if (eigs) fcg->ned = ksp->its-1; PetscFunctionReturn(0); }