12 #ifdef LIBMESH_HAVE_PETSC
13 #if !PETSC_VERSION_LESS_THAN(3, 6, 0)
14 #include "FEProblem.h"
15 #include "NonlinearSystem.h"
16 #include "libmesh/petsc_nonlinear_solver.h"
17 #include "libmesh/petsc_solver_exception.h"
32 _solver =
dynamic_cast<PetscNonlinearSolver<Real> *
>(
33 _fe_problem.getNonlinearSystem().nonlinearSolver());
36 "This line search operates only with Petsc, so Petsc must be your nonlinear solver.");
42 PetscBool changed_y = PETSC_FALSE, changed_w = PETSC_FALSE;
44 Vec X, F, Y, W, G, W1;
45 SNESLineSearch line_search;
46 PetscReal fnorm, xnorm, ynorm, gnorm;
47 PetscBool domainerror;
48 PetscReal ksp_rtol, ksp_abstol, ksp_dtol;
53 ierr = SNESGetLineSearch(snes, &line_search);
55 ierr = SNESLineSearchGetVecs(line_search, &X, &F, &Y, &W, &G);
57 ierr = SNESLineSearchGetNorms(line_search, &xnorm, &fnorm, &ynorm);
59 ierr = SNESLineSearchGetSNES(line_search, &snes);
61 ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_SUCCEEDED);
63 ierr = SNESGetKSP(snes, &ksp);
65 ierr = KSPGetTolerances(ksp, &ksp_rtol, &ksp_abstol, &ksp_dtol, &ksp_maxits);
67 ierr = VecDuplicate(W, &W1);
79 ierr = SNESLineSearchPreCheck(line_search, X, Y, &changed_y);
90 ierr = SNESComputeFunction(snes, W, F);
92 ierr = SNESGetFunctionDomainError(snes, &domainerror);
96 ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
99 ierr = VecNorm(F, NORM_2, &fnorm);
100 LIBMESH_CHKERR(ierr);
109 KSPSetTolerances(ksp,
_contact_ltol, ksp_abstol, ksp_dtol, ksp_maxits);
110 _console <<
"Contact set changed since previous non-linear iteration!\n";
113 KSPSetTolerances(ksp,
_user_ksp_rtol, ksp_abstol, ksp_dtol, ksp_maxits);
122 LIBMESH_CHKERR(ierr);
124 ierr = SNESComputeFunction(snes, W1, G);
125 LIBMESH_CHKERR(ierr);
126 ierr = SNESGetFunctionDomainError(snes, &domainerror);
127 LIBMESH_CHKERR(ierr);
130 ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
131 LIBMESH_CHKERR(ierr);
133 ierr = VecNorm(G, NORM_2, &gnorm);
134 LIBMESH_CHKERR(ierr);
138 LIBMESH_CHKERR(ierr);
140 LIBMESH_CHKERR(ierr);
152 LIBMESH_CHKERR(ierr);
154 ierr = SNESLineSearchPostCheck(line_search, X, Y, W, &changed_y, &changed_w);
155 LIBMESH_CHKERR(ierr);
159 ierr = VecWAXPY(W, -1., Y, X);
160 LIBMESH_CHKERR(ierr);
163 if (changed_w || changed_y)
165 ierr = SNESComputeFunction(snes, W, F);
166 LIBMESH_CHKERR(ierr);
167 ierr = SNESGetFunctionDomainError(snes, &domainerror);
168 LIBMESH_CHKERR(ierr);
171 ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
172 LIBMESH_CHKERR(ierr);
180 ierr = VecCopy(W, X);
181 LIBMESH_CHKERR(ierr);
183 ierr = SNESLineSearchComputeNorms(line_search);
184 LIBMESH_CHKERR(ierr);
186 ierr = VecDestroy(&W1);
187 LIBMESH_CHKERR(ierr);
192 #endif // !PETSC_VERSION_LESS_THAN(3, 3, 0)
193 #endif // LIBMESH_HAVE_PETSC