www.mooseframework.org
PetscContactLineSearch.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 
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"
18 #include <petscdm.h>
19 
21 
22 template <>
23 InputParameters
25 {
27 }
28 
29 PetscContactLineSearch::PetscContactLineSearch(const InputParameters & parameters)
30  : ContactLineSearchBase(parameters)
31 {
32  _solver = dynamic_cast<PetscNonlinearSolver<Real> *>(
33  _fe_problem.getNonlinearSystem().nonlinearSolver());
34  if (!_solver)
35  mooseError(
36  "This line search operates only with Petsc, so Petsc must be your nonlinear solver.");
37 }
38 
39 void
41 {
42  PetscBool changed_y = PETSC_FALSE, changed_w = PETSC_FALSE;
43  PetscErrorCode ierr;
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;
49  PetscInt ksp_maxits;
50  KSP ksp;
51  SNES snes = _solver->snes();
52 
53  ierr = SNESGetLineSearch(snes, &line_search);
54  LIBMESH_CHKERR(ierr);
55  ierr = SNESLineSearchGetVecs(line_search, &X, &F, &Y, &W, &G);
56  LIBMESH_CHKERR(ierr);
57  ierr = SNESLineSearchGetNorms(line_search, &xnorm, &fnorm, &ynorm);
58  LIBMESH_CHKERR(ierr);
59  ierr = SNESLineSearchGetSNES(line_search, &snes);
60  LIBMESH_CHKERR(ierr);
61  ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_SUCCEEDED);
62  LIBMESH_CHKERR(ierr);
63  ierr = SNESGetKSP(snes, &ksp);
64  LIBMESH_CHKERR(ierr);
65  ierr = KSPGetTolerances(ksp, &ksp_rtol, &ksp_abstol, &ksp_dtol, &ksp_maxits);
66  LIBMESH_CHKERR(ierr);
67  ierr = VecDuplicate(W, &W1);
68  LIBMESH_CHKERR(ierr);
69 
70  if (!_user_ksp_rtol_set)
71  {
72  _user_ksp_rtol = ksp_rtol;
73  _user_ksp_rtol_set = true;
74  }
75 
76  ++_nl_its;
77 
78  /* precheck */
79  ierr = SNESLineSearchPreCheck(line_search, X, Y, &changed_y);
80  LIBMESH_CHKERR(ierr);
81 
82  /* temporary update */
83  _contact_lambda = 1.;
84  ierr = VecWAXPY(W, -_contact_lambda, Y, X);
85  LIBMESH_CHKERR(ierr);
86 
87  /* compute residual to determine whether contact state has changed since the last non-linear
88  * residual evaluation */
89  _current_contact_state.clear();
90  ierr = SNESComputeFunction(snes, W, F);
91  LIBMESH_CHKERR(ierr);
92  ierr = SNESGetFunctionDomainError(snes, &domainerror);
93  LIBMESH_CHKERR(ierr);
94  if (domainerror)
95  {
96  ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
97  LIBMESH_CHKERR(ierr);
98  }
99  ierr = VecNorm(F, NORM_2, &fnorm);
100  LIBMESH_CHKERR(ierr);
101  std::set<dof_id_type> contact_state_stored = _current_contact_state;
102  _current_contact_state.clear();
103  printContactInfo(contact_state_stored);
104 
105  if (_affect_ltol)
106  {
107  if (contact_state_stored != _old_contact_state)
108  {
109  KSPSetTolerances(ksp, _contact_ltol, ksp_abstol, ksp_dtol, ksp_maxits);
110  _console << "Contact set changed since previous non-linear iteration!\n";
111  }
112  else
113  KSPSetTolerances(ksp, _user_ksp_rtol, ksp_abstol, ksp_dtol, ksp_maxits);
114  }
115 
116  size_t ls_its = 0;
117  while (ls_its < _allowed_lambda_cuts)
118  {
119  _contact_lambda *= 0.5;
120  /* update */
121  ierr = VecWAXPY(W1, -_contact_lambda, Y, X);
122  LIBMESH_CHKERR(ierr);
123 
124  ierr = SNESComputeFunction(snes, W1, G);
125  LIBMESH_CHKERR(ierr);
126  ierr = SNESGetFunctionDomainError(snes, &domainerror);
127  LIBMESH_CHKERR(ierr);
128  if (domainerror)
129  {
130  ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
131  LIBMESH_CHKERR(ierr);
132  }
133  ierr = VecNorm(G, NORM_2, &gnorm);
134  LIBMESH_CHKERR(ierr);
135  if (gnorm < fnorm)
136  {
137  VecCopy(G, F);
138  LIBMESH_CHKERR(ierr);
139  VecCopy(W1, W);
140  LIBMESH_CHKERR(ierr);
141  fnorm = gnorm;
142  contact_state_stored.swap(_current_contact_state);
143  _current_contact_state.clear();
144  printContactInfo(contact_state_stored);
145  ++ls_its;
146  }
147  else
148  break;
149  }
150 
151  ierr = VecScale(Y, _contact_lambda);
152  LIBMESH_CHKERR(ierr);
153  /* postcheck */
154  ierr = SNESLineSearchPostCheck(line_search, X, Y, W, &changed_y, &changed_w);
155  LIBMESH_CHKERR(ierr);
156 
157  if (changed_y)
158  {
159  ierr = VecWAXPY(W, -1., Y, X);
160  LIBMESH_CHKERR(ierr);
161  }
162 
163  if (changed_w || changed_y)
164  {
165  ierr = SNESComputeFunction(snes, W, F);
166  LIBMESH_CHKERR(ierr);
167  ierr = SNESGetFunctionDomainError(snes, &domainerror);
168  LIBMESH_CHKERR(ierr);
169  if (domainerror)
170  {
171  ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
172  LIBMESH_CHKERR(ierr);
173  }
174  contact_state_stored.swap(_current_contact_state);
175  _current_contact_state.clear();
176  printContactInfo(contact_state_stored);
177  }
178 
179  /* copy the solution over */
180  ierr = VecCopy(W, X);
181  LIBMESH_CHKERR(ierr);
182 
183  ierr = SNESLineSearchComputeNorms(line_search);
184  LIBMESH_CHKERR(ierr);
185 
186  ierr = VecDestroy(&W1);
187  LIBMESH_CHKERR(ierr);
188 
189  _old_contact_state = std::move(contact_state_stored);
190 }
191 
192 #endif // !PETSC_VERSION_LESS_THAN(3, 3, 0)
193 #endif // LIBMESH_HAVE_PETSC
validParams< ContactLineSearchBase >
InputParameters validParams< ContactLineSearchBase >()
Definition: ContactLineSearchBase.C:21
validParams< PetscContactLineSearch >
InputParameters validParams< PetscContactLineSearch >()
Definition: PetscContactLineSearch.C:24
ContactLineSearchBase::_user_ksp_rtol
Real _user_ksp_rtol
the linear tolerance set by the user in the input file
Definition: ContactLineSearchBase.h:65
ContactLineSearchBase::printContactInfo
void printContactInfo(const std::set< dof_id_type > &contact_set)
Method for printing the contact information.
Definition: ContactLineSearchBase.C:45
PetscContactLineSearch::PetscContactLineSearch
PetscContactLineSearch(const InputParameters &parameters)
Definition: PetscContactLineSearch.C:29
ContactLineSearchBase
This class implements a custom line search for use with mechanical contact.
Definition: ContactLineSearchBase.h:38
ContactLineSearchBase::_user_ksp_rtol_set
bool _user_ksp_rtol_set
Whether the user linear tolerance has been set yet in this object.
Definition: ContactLineSearchBase.h:67
PetscContactLineSearch::_solver
PetscNonlinearSolver< Real > * _solver
Definition: PetscContactLineSearch.h:44
ContactLineSearchBase::_current_contact_state
std::set< dof_id_type > _current_contact_state
The current contact set.
Definition: ContactLineSearchBase.h:60
PetscContactLineSearch.h
PetscContactLineSearch
Petsc implementation of the contact line search (based on the Petsc LineSearchShell)
Definition: PetscContactLineSearch.h:36
ContactLineSearchBase::_affect_ltol
bool _affect_ltol
Whether to modify the linear tolerance.
Definition: ContactLineSearchBase.h:79
ContactLineSearchBase::_allowed_lambda_cuts
unsigned _allowed_lambda_cuts
How many times the linsearch is allowed to cut lambda.
Definition: ContactLineSearchBase.h:73
PetscContactLineSearch::lineSearch
virtual void lineSearch() override
Definition: PetscContactLineSearch.C:40
registerMooseObject
registerMooseObject("ContactApp", PetscContactLineSearch)
ContactLineSearchBase::_contact_lambda
Real _contact_lambda
The multiplier of the newton step.
Definition: ContactLineSearchBase.h:70
ContactLineSearchBase::_contact_ltol
Real _contact_ltol
What the linear tolerance should be while the contact state is changing.
Definition: ContactLineSearchBase.h:76
ContactLineSearchBase::_old_contact_state
std::set< dof_id_type > _old_contact_state
The old contact set.
Definition: ContactLineSearchBase.h:62