www.mooseframework.org
ElementJacobianDamper.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 "ElementJacobianDamper.h"
11 #include "FEProblem.h"
12 #include "MooseMesh.h"
13 #include "DisplacedProblem.h"
14 #include "Assembly.h"
15 
16 #include "libmesh/quadrature.h" // _qrule
17 
18 registerMooseObject("TensorMechanicsApp", ElementJacobianDamper);
19 
21 
22 InputParameters
24 {
25  InputParameters params = GeneralDamper::validParams();
26  params.addClassDescription("Damper that limits the change in element Jacobians");
27  params.addParam<std::vector<VariableName>>("displacements",
28  "The nonlinear displacement variables");
29  params.addParam<Real>(
30  "max_increment",
31  0.1,
32  "The maximum permissible relative increment in the Jacobian per Newton iteration");
33  params.set<bool>("use_displaced_mesh") = true;
34  return params;
35 }
36 
37 ElementJacobianDamper::ElementJacobianDamper(const InputParameters & parameters)
38  : GeneralDamper(parameters),
39  _tid(parameters.get<THREAD_ID>("_tid")),
40  _assembly(_subproblem.assembly(_tid)),
41  _qrule(_assembly.qRule()),
42  _JxW(_assembly.JxW()),
43  _fe_problem(*parameters.getCheckedPointerParam<FEProblemBase *>("_fe_problem_base")),
44  _displaced_problem(_fe_problem.getDisplacedProblem()),
45  _max_jacobian_diff(parameters.get<Real>("max_increment"))
46 {
47  if (_displaced_problem == NULL)
48  mooseError("ElementJacobianDamper: Must use displaced problem");
49 
50  _mesh = &_displaced_problem->mesh();
51 
52  const std::vector<VariableName> & nl_vnames(getParam<std::vector<VariableName>>("displacements"));
53  _ndisp = nl_vnames.size();
54 
55  for (unsigned int i = 0; i < _ndisp; ++i)
56  {
57  _disp_var.push_back(&_sys.getFieldVariable<Real>(_tid, nl_vnames[i]));
58  _disp_incr.push_back(_disp_var.back()->increment());
59  }
60 }
61 
62 void
64 {
65 }
66 
67 Real
68 ElementJacobianDamper::computeDamping(const NumericVector<Number> & /* solution */,
69  const NumericVector<Number> & update)
70 {
71  // Maximum difference in the Jacobian for this Newton iteration
72  Real max_difference = 0.0;
73  MooseArray<Real> JxW_displaced;
74 
75  // Vector for storing the original node coordinates
76  std::vector<Point> point_copies;
77 
78  // Loop over elements in the mesh
79  for (auto & current_elem : _mesh->getMesh().active_local_element_ptr_range())
80  {
81  point_copies.clear();
82  point_copies.reserve(current_elem->n_nodes());
83 
84  // Displace nodes with current Newton increment
85  for (unsigned int i = 0; i < current_elem->n_nodes(); ++i)
86  {
87  Node & displaced_node = current_elem->node_ref(i);
88 
89  point_copies.push_back(displaced_node);
90 
91  for (unsigned int j = 0; j < _ndisp; ++j)
92  {
93  dof_id_type disp_dof_num =
94  displaced_node.dof_number(_sys.number(), _disp_var[j]->number(), 0);
95  displaced_node(j) += update(disp_dof_num);
96  }
97  }
98 
99  // Reinit element to compute Jacobian of displaced element
100  _assembly.reinit(current_elem);
101  JxW_displaced = _JxW;
102 
103  // Un-displace nodes
104  for (unsigned int i = 0; i < current_elem->n_nodes(); ++i)
105  {
106  Node & displaced_node = current_elem->node_ref(i);
107 
108  for (unsigned int j = 0; j < _ndisp; ++j)
109  displaced_node(j) = point_copies[i](j);
110  }
111 
112  // Reinit element to compute Jacobian before displacement
113  _assembly.reinit(current_elem);
114 
115  for (unsigned int qp = 0; qp < _qrule->n_points(); ++qp)
116  {
117  Real diff = std::abs(JxW_displaced[qp] - _JxW[qp]) / _JxW[qp];
118  if (diff > max_difference)
119  max_difference = diff;
120  }
121 
122  JxW_displaced.release();
123  }
124 
125  _communicator.max(max_difference);
126 
127  if (max_difference > _max_jacobian_diff)
128  return _max_jacobian_diff / max_difference;
129 
130  return 1.0;
131 }
ElementJacobianDamper
This class implements a damper that limits the change in the Jacobian of elements.
Definition: ElementJacobianDamper.h:28
ElementJacobianDamper::_displaced_problem
MooseSharedPointer< DisplacedProblem > _displaced_problem
The displaced problem.
Definition: ElementJacobianDamper.h:58
ElementJacobianDamper::_JxW
const MooseArray< Real > & _JxW
Transformed Jacobian weights.
Definition: ElementJacobianDamper.h:52
ElementJacobianDamper::_qrule
const QBase *const & _qrule
Quadrature rule.
Definition: ElementJacobianDamper.h:49
ElementJacobianDamper::computeDamping
virtual Real computeDamping(const NumericVector< Number > &, const NumericVector< Number > &update) override
Computes this Damper's damping.
Definition: ElementJacobianDamper.C:68
ElementJacobianDamper::_ndisp
unsigned int _ndisp
The number of displacement variables.
Definition: ElementJacobianDamper.h:67
registerMooseObject
registerMooseObject("TensorMechanicsApp", ElementJacobianDamper)
ElementJacobianDamper.h
ElementJacobianDamper::initialSetup
virtual void initialSetup() override
Definition: ElementJacobianDamper.C:63
defineLegacyParams
defineLegacyParams(ElementJacobianDamper)
ElementJacobianDamper::validParams
static InputParameters validParams()
Definition: ElementJacobianDamper.C:23
ElementJacobianDamper::_assembly
Assembly & _assembly
Definition: ElementJacobianDamper.h:46
validParams
InputParameters validParams()
ElementJacobianDamper::ElementJacobianDamper
ElementJacobianDamper(const InputParameters &parameters)
Definition: ElementJacobianDamper.C:37
ElementJacobianDamper::_tid
THREAD_ID _tid
Thread ID.
Definition: ElementJacobianDamper.h:45
ElementJacobianDamper::_mesh
MooseMesh * _mesh
The displaced mesh.
Definition: ElementJacobianDamper.h:61
ElementJacobianDamper::_disp_var
std::vector< MooseVariable * > _disp_var
The displacement variables.
Definition: ElementJacobianDamper.h:64
ElementJacobianDamper::_max_jacobian_diff
const Real _max_jacobian_diff
Maximum allowed relative increment in Jacobian.
Definition: ElementJacobianDamper.h:73
ElementJacobianDamper::_disp_incr
std::vector< VariableValue > _disp_incr
The current Newton increment in the displacement variables.
Definition: ElementJacobianDamper.h:70