/* KSPQCGQuadraticRoots - Computes the roots of the quadratic, ||s + step*p|| - delta = 0 such that step1 >= 0 >= step2. where delta: On entry delta must contain scalar delta. On exit delta is unchanged. step1: On entry step1 need not be specified. On exit step1 contains the non-negative root. step2: On entry step2 need not be specified. On exit step2 contains the non-positive root. C code is translated from the Fortran version of the MINPACK-2 Project, Argonne National Laboratory, Brett M. Averick and Richard G. Carter. */ static PetscErrorCode KSPQCGQuadraticRoots(Vec s,Vec p,PetscReal delta,PetscReal *step1,PetscReal *step2) { PetscReal dsq,ptp,pts,rad,sts; PetscErrorCode ierr; PetscFunctionBegin; ierr = VecDotRealPart(p,s,&pts);CHKERRQ(ierr); ierr = VecDotRealPart(p,p,&ptp);CHKERRQ(ierr); ierr = VecDotRealPart(s,s,&sts);CHKERRQ(ierr); dsq = delta*delta; rad = PetscSqrtReal((pts*pts) - ptp*(sts - dsq)); if (pts > 0.0) { *step2 = -(pts + rad)/ptp; *step1 = (sts - dsq)/(ptp * *step2); } else { *step1 = -(pts - rad)/ptp; *step2 = (sts - dsq)/(ptp * *step1); } PetscFunctionReturn(0); }
PetscErrorCode KSPSolve_QCG(KSP ksp) { /* Correpondence with documentation above: B = g = gradient, X = s = step Note: This is not coded correctly for complex arithmetic! */ KSP_QCG *pcgP = (KSP_QCG*)ksp->data; Mat Amat,Pmat; Vec W,WA,WA2,R,P,ASP,BS,X,B; PetscScalar scal,beta,rntrn,step; PetscReal q1,q2,xnorm,step1,step2,rnrm,btx,xtax; PetscReal ptasp,rtr,wtasp,bstp; PetscReal dzero = 0.0,bsnrm; PetscErrorCode ierr; PetscInt i,maxit; PC pc = ksp->pc; PCSide side; PetscBool diagonalscale; PetscFunctionBegin; ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr); if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name); if (ksp->transpose_solve) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Currently does not support transpose solve"); ksp->its = 0; maxit = ksp->max_it; WA = ksp->work[0]; R = ksp->work[1]; P = ksp->work[2]; ASP = ksp->work[3]; BS = ksp->work[4]; W = ksp->work[5]; WA2 = ksp->work[6]; X = ksp->vec_sol; B = ksp->vec_rhs; if (pcgP->delta <= dzero) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Input error: delta <= 0"); ierr = KSPGetPCSide(ksp,&side);CHKERRQ(ierr); if (side != PC_SYMMETRIC) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Requires symmetric preconditioner!"); /* Initialize variables */ ierr = VecSet(W,0.0);CHKERRQ(ierr); /* W = 0 */ ierr = VecSet(X,0.0);CHKERRQ(ierr); /* X = 0 */ ierr = PCGetOperators(pc,&Amat,&Pmat);CHKERRQ(ierr); /* Compute: BS = D^{-1} B */ ierr = PCApplySymmetricLeft(pc,B,BS);CHKERRQ(ierr); ierr = VecNorm(BS,NORM_2,&bsnrm);CHKERRQ(ierr); ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->its = 0; ksp->rnorm = bsnrm; ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); ierr = KSPLogResidualHistory(ksp,bsnrm);CHKERRQ(ierr); ierr = KSPMonitor(ksp,0,bsnrm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,0,bsnrm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) PetscFunctionReturn(0); /* Compute the initial scaled direction and scaled residual */ ierr = VecCopy(BS,R);CHKERRQ(ierr); ierr = VecScale(R,-1.0);CHKERRQ(ierr); ierr = VecCopy(R,P);CHKERRQ(ierr); ierr = VecDotRealPart(R,R,&rtr);CHKERRQ(ierr); for (i=0; i<=maxit; i++) { ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->its++; ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); /* Compute: asp = D^{-T}*A*D^{-1}*p */ ierr = PCApplySymmetricRight(pc,P,WA);CHKERRQ(ierr); ierr = KSP_MatMult(ksp,Amat,WA,WA2);CHKERRQ(ierr); ierr = PCApplySymmetricLeft(pc,WA2,ASP);CHKERRQ(ierr); /* Check for negative curvature */ ierr = VecDotRealPart(P,ASP,&ptasp);CHKERRQ(ierr); if (ptasp <= dzero) { /* Scaled negative curvature direction: Compute a step so that ||w + step*p|| = delta and QS(w + step*p) is least */ if (!i) { ierr = VecCopy(P,X);CHKERRQ(ierr); ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); scal = pcgP->delta / xnorm; ierr = VecScale(X,scal);CHKERRQ(ierr); } else { /* Compute roots of quadratic */ ierr = KSPQCGQuadraticRoots(W,P,pcgP->delta,&step1,&step2);CHKERRQ(ierr); ierr = VecDotRealPart(W,ASP,&wtasp);CHKERRQ(ierr); ierr = VecDotRealPart(BS,P,&bstp);CHKERRQ(ierr); ierr = VecCopy(W,X);CHKERRQ(ierr); q1 = step1*(bstp + wtasp + .5*step1*ptasp); q2 = step2*(bstp + wtasp + .5*step2*ptasp); if (q1 <= q2) { ierr = VecAXPY(X,step1,P);CHKERRQ(ierr); } else { ierr = VecAXPY(X,step2,P);CHKERRQ(ierr); } } pcgP->ltsnrm = pcgP->delta; /* convergence in direction of */ ksp->reason = KSP_CONVERGED_CG_NEG_CURVE; /* negative curvature */ if (!i) { ierr = PetscInfo1(ksp,"negative curvature: delta=%g\n",(double)pcgP->delta);CHKERRQ(ierr); } else { ierr = PetscInfo3(ksp,"negative curvature: step1=%g, step2=%g, delta=%g\n",(double)step1,(double)step2,(double)pcgP->delta);CHKERRQ(ierr); } } else { /* Compute step along p */ step = rtr/ptasp; ierr = VecCopy(W,X);CHKERRQ(ierr); /* x = w */ ierr = VecAXPY(X,step,P);CHKERRQ(ierr); /* x <- step*p + x */ ierr = VecNorm(X,NORM_2,&pcgP->ltsnrm);CHKERRQ(ierr); if (pcgP->ltsnrm > pcgP->delta) { /* Since the trial iterate is outside the trust region, evaluate a constrained step along p so that ||w + step*p|| = delta The positive step is always better in this case. */ if (!i) { scal = pcgP->delta / pcgP->ltsnrm; ierr = VecScale(X,scal);CHKERRQ(ierr); } else { /* Compute roots of quadratic */ ierr = KSPQCGQuadraticRoots(W,P,pcgP->delta,&step1,&step2);CHKERRQ(ierr); ierr = VecCopy(W,X);CHKERRQ(ierr); ierr = VecAXPY(X,step1,P);CHKERRQ(ierr); /* x <- step1*p + x */ } pcgP->ltsnrm = pcgP->delta; ksp->reason = KSP_CONVERGED_CG_CONSTRAINED; /* convergence along constrained step */ if (!i) { ierr = PetscInfo1(ksp,"constrained step: delta=%g\n",(double)pcgP->delta);CHKERRQ(ierr); } else { ierr = PetscInfo3(ksp,"constrained step: step1=%g, step2=%g, delta=%g\n",(double)step1,(double)step2,(double)pcgP->delta);CHKERRQ(ierr); } } else { /* Evaluate the current step */ ierr = VecCopy(X,W);CHKERRQ(ierr); /* update interior iterate */ ierr = VecAXPY(R,-step,ASP);CHKERRQ(ierr); /* r <- -step*asp + r */ ierr = VecNorm(R,NORM_2,&rnrm);CHKERRQ(ierr); ierr = PetscObjectSAWsTakeAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->rnorm = rnrm; ierr = PetscObjectSAWsGrantAccess((PetscObject)ksp);CHKERRQ(ierr); ierr = KSPLogResidualHistory(ksp,rnrm);CHKERRQ(ierr); ierr = KSPMonitor(ksp,i+1,rnrm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i+1,rnrm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) { /* convergence for */ ierr = PetscInfo3(ksp,"truncated step: step=%g, rnrm=%g, delta=%g\n",(double)PetscRealPart(step),(double)rnrm,(double)pcgP->delta);CHKERRQ(ierr); } } } if (ksp->reason) break; /* Convergence has been attained */ else { /* Compute a new AS-orthogonal direction */ ierr = VecDot(R,R,&rntrn);CHKERRQ(ierr); beta = rntrn/rtr; ierr = VecAYPX(P,beta,R);CHKERRQ(ierr); /* p <- r + beta*p */ rtr = PetscRealPart(rntrn); } } if (!ksp->reason) ksp->reason = KSP_DIVERGED_ITS; /* Unscale x */ ierr = VecCopy(X,WA2);CHKERRQ(ierr); ierr = PCApplySymmetricRight(pc,WA2,X);CHKERRQ(ierr); ierr = KSP_MatMult(ksp,Amat,X,WA);CHKERRQ(ierr); ierr = VecDotRealPart(B,X,&btx);CHKERRQ(ierr); ierr = VecDotRealPart(X,WA,&xtax);CHKERRQ(ierr); pcgP->quadratic = btx + .5*xtax; PetscFunctionReturn(0); }
PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold) { PetscErrorCode ierr; SNES_QN *qn = (SNES_QN*)snes->data; Vec W = snes->work[3]; Vec *dX = qn->U; Vec *dF = qn->V; PetscScalar *alpha = qn->alpha; PetscScalar *beta = qn->beta; PetscScalar *dXtdF = qn->dXtdF; PetscScalar *dFtdX = qn->dFtdX; PetscScalar *YtdX = qn->YtdX; /* ksp thing for jacobian scaling */ KSPConvergedReason kspreason; PetscInt k,i,j,g,lits; PetscInt m = qn->m; PetscScalar t; PetscInt l = m; PetscFunctionBegin; if (it < m) l = it; ierr = VecCopy(D,Y);CHKERRQ(ierr); if (it > 0) { k = (it - 1) % l; ierr = VecCopy(D,dF[k]);CHKERRQ(ierr); ierr = VecAXPY(dF[k], -1.0, Dold);CHKERRQ(ierr); ierr = VecCopy(X, dX[k]);CHKERRQ(ierr); ierr = VecAXPY(dX[k], -1.0, Xold);CHKERRQ(ierr); if (qn->singlereduction) { PetscScalar dFtdF; ierr = VecMDotBegin(dF[k],l,dX,dXtdF);CHKERRQ(ierr); ierr = VecMDotBegin(dX[k],l,dF,dFtdX);CHKERRQ(ierr); ierr = VecMDotBegin(Y,l,dX,YtdX);CHKERRQ(ierr); if (qn->scale_type == SNES_QN_SCALE_SHANNO) {ierr = VecDotBegin(dF[k],dF[k],&dFtdF);CHKERRQ(ierr);} ierr = VecMDotEnd(dF[k],l,dX,dXtdF);CHKERRQ(ierr); ierr = VecMDotEnd(dX[k],l,dF,dFtdX);CHKERRQ(ierr); ierr = VecMDotEnd(Y,l,dX,YtdX);CHKERRQ(ierr); if (qn->scale_type == SNES_QN_SCALE_SHANNO) { ierr = VecDotEnd(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); qn->scaling = PetscRealPart(dXtdF[k]) / PetscRealPart(dFtdF); } for (j = 0; j < l; j++) { H(k, j) = dFtdX[j]; H(j, k) = dXtdF[j]; } /* copy back over to make the computation of alpha and beta easier */ for (j = 0; j < l; j++) dXtdF[j] = H(j, j); } else { ierr = VecDot(dX[k], dF[k], &dXtdF[k]);CHKERRQ(ierr); if (qn->scale_type == SNES_QN_SCALE_SHANNO) { PetscReal dFtdF; ierr = VecDotRealPart(dF[k],dF[k],&dFtdF);CHKERRQ(ierr); qn->scaling = PetscRealPart(dXtdF[k])/dFtdF; } } if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) { ierr = SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);CHKERRQ(ierr); } } /* outward recursion starting at iteration k's update and working back */ for (i=0; i<l; i++) { k = (it-i-1)%l; if (qn->singlereduction) { /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */ t = YtdX[k]; for (j=0; j<i; j++) { g = (it-j-1)%l; t -= alpha[g]*H(k, g); } alpha[k] = t/H(k,k); } else { ierr = VecDot(dX[k],Y,&t);CHKERRQ(ierr); alpha[k] = t/dXtdF[k]; } if (qn->monitor) { ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));CHKERRQ(ierr); ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); } ierr = VecAXPY(Y,-alpha[k],dF[k]);CHKERRQ(ierr); } if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) { ierr = KSPSolve(snes->ksp,Y,W);CHKERRQ(ierr); ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr); if (kspreason < 0) { if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) { ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr); snes->reason = SNES_DIVERGED_LINEAR_SOLVE; PetscFunctionReturn(0); } } ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); snes->linear_its += lits; ierr = VecCopy(W, Y);CHKERRQ(ierr); } else { ierr = VecScale(Y, qn->scaling);CHKERRQ(ierr); } if (qn->singlereduction) { ierr = VecMDot(Y,l,dF,YtdX);CHKERRQ(ierr); } /* inward recursion starting at the first update and working forward */ for (i = 0; i < l; i++) { k = (it + i - l) % l; if (qn->singlereduction) { t = YtdX[k]; for (j = 0; j < i; j++) { g = (it + j - l) % l; t += (alpha[g] - beta[g])*H(g, k); } beta[k] = t / H(k, k); } else { ierr = VecDot(dF[k], Y, &t);CHKERRQ(ierr); beta[k] = t / dXtdF[k]; } ierr = VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);CHKERRQ(ierr); if (qn->monitor) { ierr = PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));CHKERRQ(ierr); ierr = PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);CHKERRQ(ierr); } } PetscFunctionReturn(0); }
static PetscErrorCode KSPSolve_LSQR(KSP ksp) { PetscErrorCode ierr; PetscInt i,size1,size2; PetscScalar rho,rhobar,phi,phibar,theta,c,s,tmp,tau; PetscReal beta,alpha,rnorm; Vec X,B,V,V1,U,U1,TMP,W,W2,SE,Z = NULL; Mat Amat,Pmat; MatStructure pflag; KSP_LSQR *lsqr = (KSP_LSQR*)ksp->data; PetscBool diagonalscale,nopreconditioner; PetscFunctionBegin; ierr = PCGetDiagonalScale(ksp->pc,&diagonalscale);CHKERRQ(ierr); if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name); ierr = PCGetOperators(ksp->pc,&Amat,&Pmat,&pflag);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)ksp->pc,PCNONE,&nopreconditioner);CHKERRQ(ierr); /* nopreconditioner =PETSC_FALSE; */ /* Calculate norm of right hand side */ ierr = VecNorm(ksp->vec_rhs,NORM_2,&lsqr->rhs_norm);CHKERRQ(ierr); /* mark norm of matrix with negative number to indicate it has not yet been computed */ lsqr->anorm = -1.0; /* vectors of length m, where system size is mxn */ B = ksp->vec_rhs; U = lsqr->vwork_m[0]; U1 = lsqr->vwork_m[1]; /* vectors of length n */ X = ksp->vec_sol; W = lsqr->vwork_n[0]; V = lsqr->vwork_n[1]; V1 = lsqr->vwork_n[2]; W2 = lsqr->vwork_n[3]; if (!nopreconditioner) Z = lsqr->vwork_n[4]; /* standard error vector */ SE = lsqr->se; if (SE) { ierr = VecGetSize(SE,&size1);CHKERRQ(ierr); ierr = VecGetSize(X,&size2);CHKERRQ(ierr); if (size1 != size2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Standard error vector (size %d) does not match solution vector (size %d)",size1,size2); ierr = VecSet(SE,0.0);CHKERRQ(ierr); } /* Compute initial residual, temporarily use work vector u */ if (!ksp->guess_zero) { ierr = KSP_MatMult(ksp,Amat,X,U);CHKERRQ(ierr); /* u <- b - Ax */ ierr = VecAYPX(U,-1.0,B);CHKERRQ(ierr); } else { ierr = VecCopy(B,U);CHKERRQ(ierr); /* u <- b (x is 0) */ } /* Test for nothing to do */ ierr = VecNorm(U,NORM_2,&rnorm);CHKERRQ(ierr); ierr = PetscObjectAMSTakeAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->its = 0; ksp->rnorm = rnorm; ierr = PetscObjectAMSGrantAccess((PetscObject)ksp);CHKERRQ(ierr); ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = KSPMonitor(ksp,0,rnorm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,0,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) PetscFunctionReturn(0); beta = rnorm; ierr = VecScale(U,1.0/beta);CHKERRQ(ierr); ierr = KSP_MatMultTranspose(ksp,Amat,U,V);CHKERRQ(ierr); if (nopreconditioner) { ierr = VecNorm(V,NORM_2,&alpha);CHKERRQ(ierr); } else { ierr = PCApply(ksp->pc,V,Z);CHKERRQ(ierr); ierr = VecDotRealPart(V,Z,&alpha);CHKERRQ(ierr); if (alpha <= 0.0) { ksp->reason = KSP_DIVERGED_BREAKDOWN; PetscFunctionReturn(0); } alpha = PetscSqrtReal(alpha); ierr = VecScale(Z,1.0/alpha);CHKERRQ(ierr); } ierr = VecScale(V,1.0/alpha);CHKERRQ(ierr); if (nopreconditioner) { ierr = VecCopy(V,W);CHKERRQ(ierr); } else { ierr = VecCopy(Z,W);CHKERRQ(ierr); } lsqr->arnorm = alpha * beta; phibar = beta; rhobar = alpha; i = 0; do { if (nopreconditioner) { ierr = KSP_MatMult(ksp,Amat,V,U1);CHKERRQ(ierr); } else { ierr = KSP_MatMult(ksp,Amat,Z,U1);CHKERRQ(ierr); } ierr = VecAXPY(U1,-alpha,U);CHKERRQ(ierr); ierr = VecNorm(U1,NORM_2,&beta);CHKERRQ(ierr); if (beta == 0.0) { ksp->reason = KSP_DIVERGED_BREAKDOWN; break; } ierr = VecScale(U1,1.0/beta);CHKERRQ(ierr); ierr = KSP_MatMultTranspose(ksp,Amat,U1,V1);CHKERRQ(ierr); ierr = VecAXPY(V1,-beta,V);CHKERRQ(ierr); if (nopreconditioner) { ierr = VecNorm(V1,NORM_2,&alpha);CHKERRQ(ierr); } else { ierr = PCApply(ksp->pc,V1,Z);CHKERRQ(ierr); ierr = VecDotRealPart(V1,Z,&alpha);CHKERRQ(ierr); if (alpha <= 0.0) { ksp->reason = KSP_DIVERGED_BREAKDOWN; break; } alpha = PetscSqrtReal(alpha); ierr = VecScale(Z,1.0/alpha);CHKERRQ(ierr); } ierr = VecScale(V1,1.0/alpha);CHKERRQ(ierr); rho = PetscSqrtScalar(rhobar*rhobar + beta*beta); c = rhobar / rho; s = beta / rho; theta = s * alpha; rhobar = -c * alpha; phi = c * phibar; phibar = s * phibar; tau = s * phi; ierr = VecAXPY(X,phi/rho,W);CHKERRQ(ierr); /* x <- x + (phi/rho) w */ if (SE) { ierr = VecCopy(W,W2);CHKERRQ(ierr); ierr = VecSquare(W2);CHKERRQ(ierr); ierr = VecScale(W2,1.0/(rho*rho));CHKERRQ(ierr); ierr = VecAXPY(SE, 1.0, W2);CHKERRQ(ierr); /* SE <- SE + (w^2/rho^2) */ } if (nopreconditioner) { ierr = VecAYPX(W,-theta/rho,V1);CHKERRQ(ierr); /* w <- v - (theta/rho) w */ } else { ierr = VecAYPX(W,-theta/rho,Z);CHKERRQ(ierr); /* w <- z - (theta/rho) w */ } lsqr->arnorm = alpha*PetscAbsScalar(tau); rnorm = PetscRealPart(phibar); ierr = PetscObjectAMSTakeAccess((PetscObject)ksp);CHKERRQ(ierr); ksp->its++; ksp->rnorm = rnorm; ierr = PetscObjectAMSGrantAccess((PetscObject)ksp);CHKERRQ(ierr); ierr = KSPLogResidualHistory(ksp,rnorm);CHKERRQ(ierr); ierr = KSPMonitor(ksp,i+1,rnorm);CHKERRQ(ierr); ierr = (*ksp->converged)(ksp,i+1,rnorm,&ksp->reason,ksp->cnvP);CHKERRQ(ierr); if (ksp->reason) break; SWAP(U1,U,TMP); SWAP(V1,V,TMP); i++; } while (i<ksp->max_it); if (i >= ksp->max_it && !ksp->reason) ksp->reason = KSP_DIVERGED_ITS; /* Finish off the standard error estimates */ if (SE) { tmp = 1.0; ierr = MatGetSize(Amat,&size1,&size2);CHKERRQ(ierr); if (size1 > size2) tmp = size1 - size2; tmp = rnorm / PetscSqrtScalar(tmp); ierr = VecSqrtAbs(SE);CHKERRQ(ierr); ierr = VecScale(SE,tmp);CHKERRQ(ierr); } PetscFunctionReturn(0); }
static PetscErrorCode SNESLineSearchApply_BT(SNESLineSearch linesearch) { PetscBool changed_y,changed_w; PetscErrorCode ierr; Vec X,F,Y,W,G; SNES snes; PetscReal fnorm, xnorm, ynorm, gnorm; PetscReal lambda,lambdatemp,lambdaprev,minlambda,maxstep,initslope,alpha,stol; PetscReal t1,t2,a,b,d; PetscReal f; PetscReal g,gprev; PetscBool domainerror; PetscViewer monitor; PetscInt max_its,count; SNESLineSearch_BT *bt; Mat jac; PetscErrorCode (*objective)(SNES,Vec,PetscReal*,void*); PetscFunctionBegin; ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr); ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr); ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,&maxstep,NULL,NULL,NULL,&max_its);CHKERRQ(ierr); ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr); ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr); bt = (SNESLineSearch_BT*)linesearch->data; alpha = bt->alpha; ierr = SNESGetJacobian(snes, &jac, NULL, NULL, NULL);CHKERRQ(ierr); if (!jac && !objective) SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_USER, "SNESLineSearchBT requires a Jacobian matrix"); /* precheck */ ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr); ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr); ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); if (ynorm == 0.0) { if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(monitor," Line search: Initial direction and size is 0\n");CHKERRQ(ierr); ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } ierr = VecCopy(X,W);CHKERRQ(ierr); ierr = VecCopy(F,G);CHKERRQ(ierr); ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (ynorm > maxstep) { /* Step too big, so scale back */ if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(monitor," Line search: Scaling step by %14.12e old ynorm %14.12e\n", (double)(maxstep/ynorm),(double)ynorm);CHKERRQ(ierr); ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } ierr = VecScale(Y,maxstep/(ynorm));CHKERRQ(ierr); ynorm = maxstep; } /* if the SNES has an objective set, use that instead of the function value */ if (objective) { ierr = SNESComputeObjective(snes,X,&f);CHKERRQ(ierr); } else { f = fnorm*fnorm; } /* compute the initial slope */ if (objective) { /* slope comes from the function (assumed to be the gradient of the objective */ ierr = VecDotRealPart(Y,F,&initslope);CHKERRQ(ierr); } else { /* slope comes from the normal equations */ ierr = MatMult(jac,Y,W);CHKERRQ(ierr); ierr = VecDotRealPart(F,W,&initslope);CHKERRQ(ierr); if (initslope > 0.0) initslope = -initslope; if (initslope == 0.0) initslope = -1.0; } ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); } if (snes->nfuncs >= snes->max_funcs) { ierr = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr); snes->reason = SNES_DIVERGED_FUNCTION_COUNT; ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (objective) { ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); } else { ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); if (domainerror) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (linesearch->ops->vinorm) { gnorm = fnorm; ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); } else { ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); } g = PetscSqr(gnorm); } if (PetscIsInfOrNanReal(g)) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); PetscFunctionReturn(0); } if (!objective) { ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); } if (.5*g <= .5*f + lambda*alpha*initslope) { /* Sufficient reduction or step tolerance convergence */ if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); if (!objective) { ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: fnorm %14.12e gnorm %14.12e\n", (double)fnorm, (double)gnorm);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: Using full step: obj0 %14.12e obj %14.12e\n", (double)f, (double)g);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } } else { /* Since the full step didn't work and the step is tiny, quit */ if (stol*xnorm > ynorm) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); ierr = PetscInfo2(monitor,"Aborted due to ynorm < stol*xnorm (%14.12e < %14.12e) and inadequate full step.\n",ynorm,stol*xnorm);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Fit points with quadratic */ lambdatemp = -initslope/(g - f - 2.0*lambda*initslope); lambdaprev = lambda; gprev = g; if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; if (lambdatemp <= .1*lambda) lambda = .1*lambda; else lambda = lambdatemp; ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); } if (snes->nfuncs >= snes->max_funcs) { ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr); snes->reason = SNES_DIVERGED_FUNCTION_COUNT; ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (objective) { ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); } else { ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); if (domainerror) PetscFunctionReturn(0); if (linesearch->ops->vinorm) { gnorm = fnorm; ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); } else { ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); g = gnorm*gnorm; } } if (PetscIsInfOrNanReal(g)) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); PetscFunctionReturn(0); } if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); if (!objective) { ierr = PetscViewerASCIIPrintf(monitor," Line search: gnorm after quadratic fit %14.12e\n",(double)gnorm);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: obj after quadratic fit %14.12e\n",(double)g);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } if (.5*g < .5*f + lambda*alpha*initslope) { /* sufficient reduction */ if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, lambda=%18.16e\n",(double)lambda);CHKERRQ(ierr); ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } } else { /* Fit points with cubic */ for (count = 0; count < max_its; count++) { if (lambda <= minlambda) { if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(monitor," Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr); if (!objective) { ierr = PetscViewerASCIIPrintf(monitor, " Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)fnorm, (double)gnorm, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor, " Line search: obj(0)=%18.16e, obj=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)f, (double)g, (double)ynorm, (double)minlambda, (double)lambda, (double)initslope);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { t1 = .5*(g - f) - lambda*initslope; t2 = .5*(gprev - f) - lambdaprev*initslope; a = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); b = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev); d = b*b - 3*a*initslope; if (d < 0.0) d = 0.0; if (a == 0.0) lambdatemp = -initslope/(2.0*b); else lambdatemp = (-b + PetscSqrtReal(d))/(3.0*a); } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) { lambdatemp = -initslope/(g - f - 2.0*initslope); } else SETERRQ(PetscObjectComm((PetscObject)linesearch), PETSC_ERR_SUP, "unsupported line search order for type bt"); lambdaprev = lambda; gprev = g; if (lambdatemp > .5*lambda) lambdatemp = .5*lambda; if (lambdatemp <= .1*lambda) lambda = .1*lambda; else lambda = lambdatemp; ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes,W);CHKERRQ(ierr); } if (snes->nfuncs >= snes->max_funcs) { ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr); if (!objective) { ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n", (double)fnorm,(double)gnorm,(double)ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr); } ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); snes->reason = SNES_DIVERGED_FUNCTION_COUNT; PetscFunctionReturn(0); } if (objective) { ierr = SNESComputeObjective(snes,W,&g);CHKERRQ(ierr); } else { ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); if (domainerror) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (linesearch->ops->vinorm) { gnorm = fnorm; ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); } else { ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); g = gnorm*gnorm; } } if (PetscIsInfOrNanReal(gnorm)) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); PetscFunctionReturn(0); } if (.5*g < .5*f + lambda*alpha*initslope) { /* is reduction enough? */ if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); if (!objective) { if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } else { if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratically determined step, obj %14.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } } break; } else if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); if (!objective) { if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)gnorm,(double)lambda);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } else { if (linesearch->order == SNES_LINESEARCH_ORDER_CUBIC) { ierr = PetscViewerASCIIPrintf(monitor," Line search: Cubic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: Quadratic step no good, shrinking lambda, obj %12.12e lambda=%18.16e\n",(double)g,(double)lambda);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } } } } } /* postcheck */ ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); if (changed_y) { ierr = VecWAXPY(W,-lambda,Y,X);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); } } if (changed_y || changed_w || objective) { /* recompute the function norm if the step has changed or the objective isn't the norm */ ierr = (*linesearch->ops->snesfunc)(snes,W,G);CHKERRQ(ierr); ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); if (domainerror) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } if (linesearch->ops->vinorm) { gnorm = fnorm; ierr = (*linesearch->ops->vinorm)(snes, G, W, &gnorm);CHKERRQ(ierr); } else { ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); } ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); if (PetscIsInfOrNanReal(gnorm)) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); ierr = PetscInfo(monitor,"Aborted due to Nan or Inf in function evaluation\n");CHKERRQ(ierr); PetscFunctionReturn(0); } } /* copy the solution over */ ierr = VecCopy(W, X);CHKERRQ(ierr); ierr = VecCopy(G, F);CHKERRQ(ierr); ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr); ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); ierr = SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm);CHKERRQ(ierr); PetscFunctionReturn(0); }