www.mooseframework.org
Public Member Functions | Protected Attributes | List of all members
PetscProjectSolutionOntoBounds Class Reference

Petsc implementation of the contact line search (based on the Petsc LineSearchShell) More...

#include <PetscProjectSolutionOntoBounds.h>

Inheritance diagram for PetscProjectSolutionOntoBounds:
[legend]

Public Member Functions

 PetscProjectSolutionOntoBounds (const InputParameters &parameters)
 
void initialSetup () override
 
virtual void lineSearch () override
 

Protected Attributes

NonlinearSystemBase & _nl
 
PetscNonlinearSolver< Real > * _solver
 
DisplacedProblem * _displaced_problem
 
const GeometricSearchData * _geometric_search_data
 
const std::map< std::pair< unsigned int, unsigned int >, PenetrationLocator * > * _pentration_locators
 

Detailed Description

Petsc implementation of the contact line search (based on the Petsc LineSearchShell)

Definition at line 43 of file PetscProjectSolutionOntoBounds.h.

Constructor & Destructor Documentation

◆ PetscProjectSolutionOntoBounds()

PetscProjectSolutionOntoBounds::PetscProjectSolutionOntoBounds ( const InputParameters &  parameters)

Definition at line 34 of file PetscProjectSolutionOntoBounds.C.

35  : LineSearch(parameters), _nl(_fe_problem.getNonlinearSystemBase())
36 {
37  _solver = dynamic_cast<PetscNonlinearSolver<Real> *>(
38  _fe_problem.getNonlinearSystem().nonlinearSolver());
39  if (!_solver)
40  mooseError(
41  "This line search operates only with Petsc, so Petsc must be your nonlinear solver.");
42 }

Member Function Documentation

◆ initialSetup()

void PetscProjectSolutionOntoBounds::initialSetup ( )
override

Definition at line 45 of file PetscProjectSolutionOntoBounds.C.

46 {
47  _displaced_problem = _fe_problem.getDisplacedProblem().get();
48  if (!_displaced_problem)
49  mooseError("This line search only makes sense in a displaced context");
50 }

◆ lineSearch()

void PetscProjectSolutionOntoBounds::lineSearch ( )
overridevirtual

Definition at line 53 of file PetscProjectSolutionOntoBounds.C.

54 {
55  PetscBool changed_y = PETSC_FALSE;
56  PetscErrorCode ierr;
57  Vec X, F, Y, W, G;
58  SNESLineSearch line_search;
59  PetscReal fnorm, xnorm, ynorm;
60  PetscBool domainerror;
61  SNES snes = _solver->snes();
62 
63  ierr = SNESGetLineSearch(snes, &line_search);
64  LIBMESH_CHKERR(ierr);
65  ierr = SNESLineSearchGetVecs(line_search, &X, &F, &Y, &W, &G);
66  LIBMESH_CHKERR(ierr);
67  ierr = SNESLineSearchGetNorms(line_search, &xnorm, &fnorm, &ynorm);
68  LIBMESH_CHKERR(ierr);
69  ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_SUCCEEDED);
70  LIBMESH_CHKERR(ierr);
71 
72  ierr = SNESLineSearchPreCheck(line_search, X, Y, &changed_y);
73  LIBMESH_CHKERR(ierr);
74 
75  // apply the full newton step
76  ierr = VecWAXPY(W, -1., Y, X);
77  LIBMESH_CHKERR(ierr);
78 
79  {
80  PetscVector<Number> solution(W, this->comm());
81 
82  _nl.setSolution(solution);
83 
84  // Displace the mesh and update the displaced geometric search objects
85  _displaced_problem->updateMesh();
86  }
87 
88  // A reference to the displaced MeshBase object
89  const auto & mesh = _displaced_problem->mesh().getMesh();
90 
91  const auto & pen_locs = _displaced_problem->geomSearchData().getPenetrationLocators();
92 
93  // Keep track of the slave nodes that we push back onto the master face. We'll eventually check to
94  // make sure that we didn't have a corner node at the intersection of two slave faces that we
95  // tried to displace twice. As this stands now this won't cover the case wherethe intersection
96  // happens only across processes
97  std::set<dof_id_type> nodes_displaced;
98 
99  std::vector<unsigned int> disp_nums;
100 
101  // Generate the displaced variable numbers
102  for (const auto & disp_name : _displaced_problem->getDisplacementVarNames())
103  disp_nums.push_back(_nl.system().variable_number(disp_name));
104 
105  for (const auto & pen_loc : pen_locs)
106  for (const auto & pinfo_pair : pen_loc.second->_penetration_info)
107  {
108  auto node_id = pinfo_pair.first;
109  auto pen_info = pinfo_pair.second;
110 
111  // We have penetration
112  if (pen_info->_distance > 0)
113  {
114 // Avoid warning in optimized modes about unused variables
115 #ifndef NDEBUG
116  // Make sure we haven't done this node before
117  auto pair = nodes_displaced.insert(node_id);
118  mooseAssert(pair.second, "Node id " << node_id << " has already been displaced");
119 #endif
120 
121  const auto & node = mesh.node_ref(node_id);
122 
123  // If this is not a local node, we will let displacement happen on another process
124  if (node.processor_id() != this->processor_id())
125  continue;
126 
127  // The vector that we need to displace by
128  auto required_solution_change = pen_info->_distance * pen_info->_normal;
129 
130  unsigned component = 0;
131  std::vector<PetscInt> indices;
132  std::vector<PetscScalar> values;
133  for (auto disp_num : disp_nums)
134  {
135  auto dof_number = node.dof_number(/*sys=*/0, disp_num, /*component=*/0);
136  indices.push_back(static_cast<PetscInt>(dof_number));
137  values.push_back(static_cast<PetscScalar>(required_solution_change(component++)));
138  }
139  ierr = VecSetValues(
140  W, static_cast<PetscInt>(indices.size()), indices.data(), values.data(), ADD_VALUES);
141  LIBMESH_CHKERR(ierr);
142  }
143  }
144 
145  ierr = VecAssemblyBegin(W);
146  LIBMESH_CHKERR(ierr);
147  ierr = VecAssemblyEnd(W);
148  LIBMESH_CHKERR(ierr);
149 
150  ierr = SNESComputeFunction(snes, W, F);
151  LIBMESH_CHKERR(ierr);
152  ierr = SNESGetFunctionDomainError(snes, &domainerror);
153  LIBMESH_CHKERR(ierr);
154  if (domainerror)
155  {
156  ierr = SNESLineSearchSetReason(line_search, SNES_LINESEARCH_FAILED_DOMAIN);
157  LIBMESH_CHKERR(ierr);
158  }
159  ierr = VecNorm(F, NORM_2, &fnorm);
160  LIBMESH_CHKERR(ierr);
161 
162  ierr = VecCopy(W, X);
163  LIBMESH_CHKERR(ierr);
164 
165  ierr = SNESLineSearchComputeNorms(line_search);
166  LIBMESH_CHKERR(ierr);
167 }

Member Data Documentation

◆ _displaced_problem

DisplacedProblem* PetscProjectSolutionOntoBounds::_displaced_problem
protected

Definition at line 54 of file PetscProjectSolutionOntoBounds.h.

Referenced by initialSetup(), and lineSearch().

◆ _geometric_search_data

const GeometricSearchData* PetscProjectSolutionOntoBounds::_geometric_search_data
protected

Definition at line 55 of file PetscProjectSolutionOntoBounds.h.

◆ _nl

NonlinearSystemBase& PetscProjectSolutionOntoBounds::_nl
protected

Definition at line 52 of file PetscProjectSolutionOntoBounds.h.

Referenced by lineSearch().

◆ _pentration_locators

const std::map<std::pair<unsigned int, unsigned int>, PenetrationLocator *>* PetscProjectSolutionOntoBounds::_pentration_locators
protected

Definition at line 57 of file PetscProjectSolutionOntoBounds.h.

◆ _solver

PetscNonlinearSolver<Real>* PetscProjectSolutionOntoBounds::_solver
protected

Definition at line 53 of file PetscProjectSolutionOntoBounds.h.

Referenced by lineSearch(), and PetscProjectSolutionOntoBounds().


The documentation for this class was generated from the following files:
PetscProjectSolutionOntoBounds::_solver
PetscNonlinearSolver< Real > * _solver
Definition: PetscProjectSolutionOntoBounds.h:53
PetscProjectSolutionOntoBounds::_displaced_problem
DisplacedProblem * _displaced_problem
Definition: PetscProjectSolutionOntoBounds.h:54
PetscProjectSolutionOntoBounds::_nl
NonlinearSystemBase & _nl
Definition: PetscProjectSolutionOntoBounds.h:52
MaterialTensorCalculatorTools::component
Real component(const SymmTensor &symm_tensor, unsigned int index)
Definition: MaterialTensorCalculatorTools.C:16