https://mooseframework.inl.gov
PetscContactLineSearch.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "PetscContactLineSearch.h"
11 #include "FEProblem.h"
12 #include "NonlinearSystem.h"
13 #include "libmesh/petsc_nonlinear_solver.h"
14 #include "libmesh/petsc_solver_exception.h"
15 #include <petscdm.h>
16 
18 
21 {
23 }
24 
26  : ContactLineSearchBase(parameters)
27 {
28  _solver = dynamic_cast<PetscNonlinearSolver<Real> *>(
29  _fe_problem.getNonlinearSystem(/*nl_sys_num=*/0).nonlinearSolver());
30  if (!_solver)
31  mooseError(
32  "This line search operates only with Petsc, so Petsc must be your nonlinear solver.");
33 }
34 
35 void
37 {
38  PetscBool changed_y = PETSC_FALSE, changed_w = PETSC_FALSE;
39  Vec X, F, Y, W, G, W1;
40  SNESLineSearch line_search;
41  PetscReal fnorm, xnorm, ynorm, gnorm;
42  PetscBool domainerror;
43  PetscReal ksp_rtol, ksp_abstol, ksp_dtol;
44  PetscInt ksp_maxits;
45  KSP ksp;
46  SNES snes = _solver->snes();
47 
48  LibmeshPetscCall(SNESGetLineSearch(snes, &line_search));
49  LibmeshPetscCall(SNESLineSearchGetVecs(line_search, &X, &F, &Y, &W, &G));
50  LibmeshPetscCall(SNESLineSearchGetNorms(line_search, &xnorm, &fnorm, &ynorm));
51  LibmeshPetscCall(SNESLineSearchGetSNES(line_search, &snes));
52  LibmeshPetscCall(SNESLineSearchSetReason(line_search, SNES_LINESEARCH_SUCCEEDED));
53  LibmeshPetscCall(SNESGetKSP(snes, &ksp));
54  LibmeshPetscCall(KSPGetTolerances(ksp, &ksp_rtol, &ksp_abstol, &ksp_dtol, &ksp_maxits));
55  LibmeshPetscCall(VecDuplicate(W, &W1));
56 
57  if (!_user_ksp_rtol_set)
58  {
59  _user_ksp_rtol = ksp_rtol;
60  _user_ksp_rtol_set = true;
61  }
62 
63  ++_nl_its;
64 
65  /* precheck */
66  LibmeshPetscCall(SNESLineSearchPreCheck(line_search, X, Y, &changed_y));
67 
68  /* temporary update */
69  _contact_lambda = 1.;
70  LibmeshPetscCall(VecWAXPY(W, -_contact_lambda, Y, X));
71 
72  /* compute residual to determine whether contact state has changed since the last non-linear
73  * residual evaluation */
74  _current_contact_state.clear();
75  LibmeshPetscCall(SNESComputeFunction(snes, W, F));
76  LibmeshPetscCall(SNESGetFunctionDomainError(snes, &domainerror));
77  if (domainerror)
78  LibmeshPetscCall(SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN));
79 
80  LibmeshPetscCall(VecNorm(F, NORM_2, &fnorm));
81  std::set<dof_id_type> contact_state_stored = _current_contact_state;
82  _current_contact_state.clear();
83  printContactInfo(contact_state_stored);
84 
85  if (_affect_ltol)
86  {
87  if (contact_state_stored != _old_contact_state)
88  {
89  LibmeshPetscCall(KSPSetTolerances(ksp, _contact_ltol, ksp_abstol, ksp_dtol, ksp_maxits));
90  _console << "Contact set changed since previous non-linear iteration!" << std::endl;
91  }
92  else
93  LibmeshPetscCall(KSPSetTolerances(ksp, _user_ksp_rtol, ksp_abstol, ksp_dtol, ksp_maxits));
94  }
95 
96  size_t ls_its = 0;
97  while (ls_its < _allowed_lambda_cuts)
98  {
99  _contact_lambda *= 0.5;
100  /* update */
101  LibmeshPetscCall(VecWAXPY(W1, -_contact_lambda, Y, X));
102 
103  LibmeshPetscCall(SNESComputeFunction(snes, W1, G));
104  LibmeshPetscCall(SNESGetFunctionDomainError(snes, &domainerror));
105  if (domainerror)
106  LibmeshPetscCall(SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN));
107 
108  LibmeshPetscCall(VecNorm(G, NORM_2, &gnorm));
109  if (gnorm < fnorm)
110  {
111  LibmeshPetscCall(VecCopy(G, F));
112  LibmeshPetscCall(VecCopy(W1, W));
113  fnorm = gnorm;
114  contact_state_stored.swap(_current_contact_state);
115  _current_contact_state.clear();
116  printContactInfo(contact_state_stored);
117  ++ls_its;
118  }
119  else
120  break;
121  }
122 
123  LibmeshPetscCall(VecScale(Y, _contact_lambda));
124  /* postcheck */
125  LibmeshPetscCall(SNESLineSearchPostCheck(line_search, X, Y, W, &changed_y, &changed_w));
126 
127  if (changed_y)
128  LibmeshPetscCall(VecWAXPY(W, -1., Y, X));
129 
130  if (changed_w || changed_y)
131  {
132  LibmeshPetscCall(SNESComputeFunction(snes, W, F));
133  LibmeshPetscCall(SNESGetFunctionDomainError(snes, &domainerror));
134  if (domainerror)
135  LibmeshPetscCall(SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN));
136 
137  contact_state_stored.swap(_current_contact_state);
138  _current_contact_state.clear();
139  printContactInfo(contact_state_stored);
140  }
141 
142  /* copy the solution over */
143  LibmeshPetscCall(VecCopy(W, X));
144 
145  LibmeshPetscCall(SNESLineSearchComputeNorms(line_search));
146 
147  LibmeshPetscCall(VecDestroy(&W1));
148 
149  _old_contact_state = std::move(contact_state_stored);
150 }
bool _affect_ltol
Whether to modify the linear tolerance.
SNES snes(const char *name=nullptr)
virtual libMesh::NonlinearSolver< Number > * nonlinearSolver() override
std::set< dof_id_type > _old_contact_state
The old contact set.
Real _contact_ltol
What the linear tolerance should be while the contact state is changing.
std::set< dof_id_type > _current_contact_state
The current contact set.
virtual void lineSearch() override
FEProblem & _fe_problem
static const std::string F
Definition: NS.h:165
static const std::string G
Definition: NS.h:166
registerMooseObject("ContactApp", PetscContactLineSearch)
This class implements a custom line search for use with mechanical contact.
Petsc implementation of the contact line search (based on the Petsc LineSearchShell) ...
Real _contact_lambda
The multiplier of the newton step.
static InputParameters validParams()
PetscContactLineSearch(const InputParameters &parameters)
unsigned _allowed_lambda_cuts
How many times the linsearch is allowed to cut lambda.
void printContactInfo(const std::set< dof_id_type > &contact_set)
Method for printing the contact information.
void mooseError(Args &&... args) const
const ConsoleStream _console
Real _user_ksp_rtol
the linear tolerance set by the user in the input file
PetscNonlinearSolver< Real > * _solver
size_t _nl_its
static InputParameters validParams()
virtual NonlinearSystem & getNonlinearSystem(const unsigned int nl_sys_num) override
bool _user_ksp_rtol_set
Whether the user linear tolerance has been set yet in this object.