https://mooseframework.inl.gov
ReferenceElementJacobianDamper.h
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 #pragma once
11 
12 // Moose Includes
13 #include "GeneralDamper.h"
14 #include "MooseVariable.h"
15 #include "RankTwoTensorForward.h"
16 #include "Coupleable.h"
17 #include "BlockRestrictable.h"
18 
19 // Forward Declarations
20 class FEProblemBase;
21 class MooseMesh;
22 
28  public Coupleable,
29  public BlockRestrictable
30 {
31 public:
33 
35 
36  virtual void initialSetup() override {}
37 
38  virtual Real computeDamping(const NumericVector<Number> & solution,
39  const NumericVector<Number> & update) override;
40 
41 protected:
44 
47 
50 
53 
55  const QBase * const & _qrule;
56 
58  const unsigned int _ndisp;
59 
61  std::vector<unsigned int> _disp_num;
62 
64  std::vector<const VariablePhiGradient *> _grad_phi;
65 
66 private:
68  void computeGradDisp(const Elem * elem,
69  const NumericVector<Number> & solution,
70  const NumericVector<Number> & update);
71 
73  std::vector<std::vector<RealVectorValue>> _grad_disp;
74 
76  std::vector<std::vector<RealVectorValue>> _grad_disp_update;
77 };
void computeGradDisp(const Elem *elem, const NumericVector< Number > &solution, const NumericVector< Number > &update)
Fill the displacement gradients.
MooseMesh & _mesh
The undisplaced mesh.
Assembly & _assembly
The undisplaced assembly.
virtual Real computeDamping(const NumericVector< Number > &solution, const NumericVector< Number > &update) override
const Real _max_jacobian_diff
Maximum allowed relative increment in Jacobian.
This class implements a damper that limits the change in the Jacobian of elements without relying on ...
std::vector< unsigned int > _disp_num
The displacement variable numbers.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
ReferenceElementJacobianDamper(const InputParameters &parameters)
const InputParameters & parameters() const
std::vector< std::vector< RealVectorValue > > _grad_disp
The current displacement gradients.
std::vector< const VariablePhiGradient * > _grad_phi
shape function gradients
const unsigned int _ndisp
Number of displacement variables.
const QBase *const & _qrule
Quadrature rule.
unsigned int THREAD_ID
std::vector< std::vector< RealVectorValue > > _grad_disp_update
The displacement gradients after this update.