/*@C SNESUpdateCheckJacobian - Checks each Jacobian computed by the nonlinear solver comparing the users function with a finite difference computation. Options Database: + -snes_check_jacobian - use this every time SNESSolve() is called - -snes_check_jacobian_view - Display difference between approximate and hand-coded Jacobian Level: intermediate .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESSolve() @*/ PetscErrorCode SNESUpdateCheckJacobian(SNES snes,PetscInt it) { Mat A = snes->jacobian,B; Vec x = snes->vec_sol,f = snes->vec_func,f1 = snes->vec_sol_update; PetscErrorCode ierr; PetscReal nrm,gnorm; PetscErrorCode (*objective)(SNES,Vec,PetscReal*,void*); void *ctx; PetscReal fnorm,f1norm,dnorm; PetscInt m,n,M,N; PetscBool complete_print = PETSC_FALSE; void *functx; PetscViewer viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); PetscFunctionBegin; ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_check_jacobian_view",&complete_print);CHKERRQ(ierr); if (A != snes->jacobian_pre) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot check Jacobian with alternative preconditioner"); ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer," Testing hand-coded Jacobian, if the ratio is O(1.e-8), the hand-coded Jacobian is probably correct.\n");CHKERRQ(ierr); if (!complete_print) { ierr = PetscViewerASCIIPrintf(viewer," Run with -snes_check_jacobian_view [viewer][:filename][:format] to show difference of hand-coded and finite difference Jacobian.\n");CHKERRQ(ierr); } /* compute both versions of Jacobian */ ierr = SNESComputeJacobian(snes,x,A,A);CHKERRQ(ierr); ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,M,N);CHKERRQ(ierr); ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); ierr = MatSetUp(B);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); ierr = SNESGetFunction(snes,NULL,NULL,&functx);CHKERRQ(ierr); ierr = SNESComputeJacobianDefault(snes,x,B,B,functx);CHKERRQ(ierr); if (complete_print) { ierr = PetscViewerASCIIPrintf(viewer," Finite difference Jacobian\n");CHKERRQ(ierr); ierr = MatViewFromOptions(B,((PetscObject)snes)->prefix,"-snes_check_jacobian_view");CHKERRQ(ierr); } /* compare */ ierr = MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); ierr = MatNorm(B,NORM_FROBENIUS,&nrm);CHKERRQ(ierr); ierr = MatNorm(A,NORM_FROBENIUS,&gnorm);CHKERRQ(ierr); if (complete_print) { ierr = PetscViewerASCIIPrintf(viewer," Hand-coded Jacobian\n");CHKERRQ(ierr); ierr = MatViewFromOptions(A,((PetscObject)snes)->prefix,"-snes_check_jacobian_view");CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer," Hand-coded minus finite difference Jacobian\n");CHKERRQ(ierr); ierr = MatViewFromOptions(B,((PetscObject)snes)->prefix,"-snes_check_jacobian_view");CHKERRQ(ierr); } if (!gnorm) gnorm = 1; /* just in case */ ierr = PetscViewerASCIIPrintf(viewer," %g = ||J - Jfd||//J|| %g = ||J - Jfd||\n",(double)(nrm/gnorm),(double)nrm);CHKERRQ(ierr); ierr = SNESGetObjective(snes,&objective,&ctx);CHKERRQ(ierr); if (objective) { ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); ierr = VecNorm(f,NORM_2,&fnorm);CHKERRQ(ierr); if (complete_print) { ierr = PetscViewerASCIIPrintf(viewer," Hand-coded Objective Function \n");CHKERRQ(ierr); ierr = VecView(f,viewer);CHKERRQ(ierr); } ierr = SNESObjectiveComputeFunctionDefaultFD(snes,x,f1,NULL);CHKERRQ(ierr); ierr = VecNorm(f1,NORM_2,&f1norm);CHKERRQ(ierr); if (complete_print) { ierr = PetscViewerASCIIPrintf(viewer," Finite-Difference Objective Function\n");CHKERRQ(ierr); ierr = VecView(f1,viewer);CHKERRQ(ierr); } /* compare the two */ ierr = VecAXPY(f,-1.0,f1);CHKERRQ(ierr); ierr = VecNorm(f,NORM_2,&dnorm);CHKERRQ(ierr); if (!fnorm) fnorm = 1.; ierr = PetscViewerASCIIPrintf(viewer," %g = Norm of objective function ratio %g = difference\n",dnorm/fnorm,dnorm);CHKERRQ(ierr); if (complete_print) { ierr = PetscViewerASCIIPrintf(viewer," Difference\n");CHKERRQ(ierr); ierr = VecView(f,viewer);CHKERRQ(ierr); } } ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)snes)->tablevel);CHKERRQ(ierr); ierr = MatDestroy(&B);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); }
PetscErrorCode SNESSolve_Test(SNES snes) { Mat A = snes->jacobian,B; Vec x = snes->vec_sol,f = snes->vec_func,f1 = snes->vec_sol_update; PetscErrorCode ierr; PetscInt i; PetscReal nrm,gnorm; SNES_Test *neP = (SNES_Test*)snes->data; PetscErrorCode (*objective)(SNES,Vec,PetscReal*,void*); void *ctx; PetscReal fnorm,f1norm,dnorm; PetscFunctionBegin; if (A != snes->jacobian_pre) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot test with alternative preconditioner"); ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Testing hand-coded Jacobian, if the ratio is\n");CHKERRQ(ierr); ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"O(1.e-8), the hand-coded Jacobian is probably correct.\n");CHKERRQ(ierr); if (!neP->complete_print) { ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Run with -snes_test_display to show difference\n");CHKERRQ(ierr); ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"of hand-coded and finite difference Jacobian.\n");CHKERRQ(ierr); } for (i=0; i<3; i++) { void *functx; static const char *const loc[] = {"user-defined state","constant state -1.0","constant state 1.0"}; PetscInt m,n,M,N; if (i == 1) { ierr = VecSet(x,-1.0);CHKERRQ(ierr); } else if (i == 2) { ierr = VecSet(x,1.0);CHKERRQ(ierr); } /* evaluate the function at this point because SNESComputeJacobianDefaultColor() assumes that the function has been evaluated and put into snes->vec_func */ ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); if (snes->domainerror) { ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Domain error at %s\n",loc[i]);CHKERRQ(ierr); snes->domainerror = PETSC_FALSE; continue; } /* compute both versions of Jacobian */ ierr = SNESComputeJacobian(snes,x,A,A);CHKERRQ(ierr); ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); ierr = MatSetSizes(B,m,n,M,N);CHKERRQ(ierr); ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); ierr = MatSetUp(B);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); ierr = SNESGetFunction(snes,NULL,NULL,&functx);CHKERRQ(ierr); ierr = SNESComputeJacobianDefault(snes,x,B,B,functx);CHKERRQ(ierr); if (neP->complete_print) { MPI_Comm comm; PetscViewer viewer; ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Finite difference Jacobian (%s)\n",loc[i]);CHKERRQ(ierr); ierr = PetscObjectGetComm((PetscObject)B,&comm);CHKERRQ(ierr); ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); ierr = MatView(B,viewer);CHKERRQ(ierr); } /* compare */ ierr = MatAYPX(B,-1.0,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); ierr = MatNorm(B,NORM_FROBENIUS,&nrm);CHKERRQ(ierr); ierr = MatNorm(A,NORM_FROBENIUS,&gnorm);CHKERRQ(ierr); if (neP->complete_print) { MPI_Comm comm; PetscViewer viewer; ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Hand-coded Jacobian (%s)\n",loc[i]);CHKERRQ(ierr); ierr = PetscObjectGetComm((PetscObject)B,&comm);CHKERRQ(ierr); ierr = PetscViewerASCIIGetStdout(comm,&viewer);CHKERRQ(ierr); ierr = MatView(A,viewer);CHKERRQ(ierr); ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Hand-coded minus finite difference Jacobian (%s)\n",loc[i]);CHKERRQ(ierr); ierr = MatView(B,viewer);CHKERRQ(ierr); } if (!gnorm) gnorm = 1; /* just in case */ ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Norm of matrix ratio %g difference %g (%s)\n",(double)(nrm/gnorm),(double)nrm,loc[i]);CHKERRQ(ierr); ierr = SNESGetObjective(snes,&objective,&ctx);CHKERRQ(ierr); if (objective) { ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr); ierr = VecNorm(f,NORM_2,&fnorm);CHKERRQ(ierr); if (neP->complete_print) { PetscViewer viewer; ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Hand-coded Function (%s)\n",loc[i]);CHKERRQ(ierr); ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);CHKERRQ(ierr); ierr = VecView(f,viewer);CHKERRQ(ierr); } ierr = SNESObjectiveComputeFunctionDefaultFD(snes,x,f1,NULL);CHKERRQ(ierr); ierr = VecNorm(f1,NORM_2,&f1norm);CHKERRQ(ierr); if (neP->complete_print) { PetscViewer viewer; ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Finite-Difference Function (%s)\n",loc[i]);CHKERRQ(ierr); ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);CHKERRQ(ierr); ierr = VecView(f1,viewer);CHKERRQ(ierr); } /* compare the two */ ierr = VecAXPY(f,-1.0,f1);CHKERRQ(ierr); ierr = VecNorm(f,NORM_2,&dnorm);CHKERRQ(ierr); if (!fnorm) fnorm = 1.; ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Norm of function ratio %g difference %g (%s)\n",dnorm/fnorm,dnorm,loc[i]);CHKERRQ(ierr); if (neP->complete_print) { PetscViewer viewer; ierr = PetscPrintf(PetscObjectComm((PetscObject)snes),"Difference (%s)\n",loc[i]);CHKERRQ(ierr); ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)snes),&viewer);CHKERRQ(ierr); ierr = VecView(f,viewer);CHKERRQ(ierr); } } ierr = MatDestroy(&B);CHKERRQ(ierr); } /* Abort after the first iteration due to the jacobian not being valid. */ SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESTest aborts after Jacobian test"); PetscFunctionReturn(0); }
static PetscErrorCode SNESLineSearchApply_L2(SNESLineSearch linesearch) { PetscBool changed_y, changed_w; PetscErrorCode ierr; Vec X; Vec F; Vec Y; Vec W; SNES snes; PetscReal gnorm; PetscReal ynorm; PetscReal xnorm; PetscReal steptol, maxstep, rtol, atol, ltol; PetscViewer monitor; PetscBool domainerror; PetscReal lambda, lambda_old, lambda_mid, lambda_update, delLambda; PetscReal fnrm, fnrm_old, fnrm_mid; PetscReal delFnrm, delFnrm_old, del2Fnrm; PetscInt i, max_its; PetscErrorCode (*objective)(SNES,Vec,PetscReal*,void*); PetscFunctionBegin; ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL);CHKERRQ(ierr); ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr); ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr); ierr = SNESLineSearchSetSuccess(linesearch, PETSC_TRUE);CHKERRQ(ierr); ierr = SNESLineSearchGetTolerances(linesearch, &steptol, &maxstep, &rtol, &atol, <ol, &max_its);CHKERRQ(ierr); ierr = SNESLineSearchGetMonitor(linesearch, &monitor);CHKERRQ(ierr); ierr = SNESGetObjective(snes,&objective,NULL);CHKERRQ(ierr); /* precheck */ ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr); lambda_old = 0.0; if (!objective) { fnrm_old = gnorm*gnorm; } else { ierr = SNESComputeObjective(snes,X,&fnrm_old);CHKERRQ(ierr); } lambda_mid = 0.5*(lambda + lambda_old); for (i = 0; i < max_its; i++) { /* compute the norm at the midpoint */ ierr = VecCopy(X, W);CHKERRQ(ierr); ierr = VecAXPY(W, -lambda_mid, Y);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); } if (!objective) { ierr = (*linesearch->ops->snesfunc)(snes, W, F);CHKERRQ(ierr); if (linesearch->ops->vinorm) { fnrm_mid = gnorm; ierr = (*linesearch->ops->vinorm)(snes, F, W, &fnrm_mid);CHKERRQ(ierr); } else { ierr = VecNorm(F,NORM_2,&fnrm_mid);CHKERRQ(ierr); } /* compute the norm at lambda */ ierr = VecCopy(X, W);CHKERRQ(ierr); ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); } ierr = (*linesearch->ops->snesfunc)(snes, W, F);CHKERRQ(ierr); if (linesearch->ops->vinorm) { fnrm = gnorm; ierr = (*linesearch->ops->vinorm)(snes, F, W, &fnrm);CHKERRQ(ierr); } else { ierr = VecNorm(F,NORM_2,&fnrm);CHKERRQ(ierr); } fnrm_mid = fnrm_mid*fnrm_mid; fnrm = fnrm*fnrm; } else { /* compute the objective at the midpoint */ ierr = VecCopy(X, W);CHKERRQ(ierr); ierr = VecAXPY(W, -lambda_mid, Y);CHKERRQ(ierr); ierr = SNESComputeObjective(snes,W,&fnrm_mid);CHKERRQ(ierr); /* compute the objective at the midpoint */ ierr = VecCopy(X, W);CHKERRQ(ierr); ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); ierr = SNESComputeObjective(snes,W,&fnrm);CHKERRQ(ierr); } /* this gives us the derivatives at the endpoints -- compute them from here a = x - a0 p_0(x) = (x / dA - 1)(2x / dA - 1) p_1(x) = 4(x / dA)(1 - x / dA) p_2(x) = (x / dA)(2x / dA - 1) dp_0[0] / dx = 3 / dA dp_1[0] / dx = -4 / dA dp_2[0] / dx = 1 / dA dp_0[dA] / dx = -1 / dA dp_1[dA] / dx = 4 / dA dp_2[dA] / dx = -3 / dA d^2p_0[0] / dx^2 = 4 / dA^2 d^2p_1[0] / dx^2 = -8 / dA^2 d^2p_2[0] / dx^2 = 4 / dA^2 */ delLambda = lambda - lambda_old; delFnrm = (3.*fnrm - 4.*fnrm_mid + 1.*fnrm_old) / delLambda; delFnrm_old = (-3.*fnrm_old + 4.*fnrm_mid -1.*fnrm) / delLambda; del2Fnrm = (delFnrm - delFnrm_old) / delLambda; if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); if (!objective) { ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g, %g], fnorms = [%g, %g, %g]\n", (double)lambda, (double)lambda_mid, (double)lambda_old, (double)PetscSqrtReal(fnrm), (double)PetscSqrtReal(fnrm_mid), (double)PetscSqrtReal(fnrm_old));CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(monitor," Line search: lambdas = [%g, %g, %g], obj = [%g, %g, %g]\n", (double)lambda, (double)lambda_mid, (double)lambda_old, (double)fnrm, (double)fnrm_mid, (double)fnrm_old);CHKERRQ(ierr); } ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } /* compute the search direction -- always go downhill */ if (del2Fnrm > 0.) lambda_update = lambda - delFnrm / del2Fnrm; else lambda_update = lambda + delFnrm / del2Fnrm; if (lambda_update < steptol) lambda_update = 0.5*(lambda + lambda_old); if (PetscIsInfOrNanScalar(lambda_update)) break; if (lambda_update > maxstep) break; /* compute the new state of the line search */ lambda_old = lambda; lambda = lambda_update; fnrm_old = fnrm; lambda_mid = 0.5*(lambda + lambda_old); } /* construct the solution */ ierr = VecCopy(X, W);CHKERRQ(ierr); ierr = VecAXPY(W, -lambda, Y);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, W);CHKERRQ(ierr); } /* postcheck */ ierr = SNESLineSearchPostCheck(linesearch,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); if (changed_y) { ierr = VecAXPY(X, -lambda, Y);CHKERRQ(ierr); if (linesearch->ops->viproject) { ierr = (*linesearch->ops->viproject)(snes, X);CHKERRQ(ierr); } } else { ierr = VecCopy(W, X);CHKERRQ(ierr); } ierr = (*linesearch->ops->snesfunc)(snes,X,F);CHKERRQ(ierr); ierr = SNESGetFunctionDomainError(snes, &domainerror);CHKERRQ(ierr); if (domainerror) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr); ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm);CHKERRQ(ierr); if (monitor) { ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(monitor," Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm);CHKERRQ(ierr); ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr); } if (lambda <= steptol) { ierr = SNESLineSearchSetSuccess(linesearch, PETSC_FALSE);CHKERRQ(ierr); } PetscFunctionReturn(0); }