https://mooseframework.inl.gov
CZMInterfaceKernelTotalLagrangian.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 
11 
13 
16 {
18  params.addClassDescription("Calculate residual contribution for balancing the traction across an "
19  "interface (used in the cohesive zone method).");
20  params.set<std::string>("traction_global_name") = "PK1traction";
21  return params;
22 }
23 
25  const InputParameters & parameters)
26  : CZMInterfaceKernelBase(parameters),
27  _dPK1traction_dF(getMaterialPropertyByName<RankThreeTensor>(_base_name + "dPK1traction_dF"))
28 {
29 }
30 
31 Real
33  const unsigned int & component_j, const Moose::DGJacobianType & type) const
34 {
35  Real jacsd = _dtraction_djump_global[_qp](_component, component_j);
36  Real jac = 0;
37 
38  switch (type)
39  {
40  case Moose::ElementElement: // Residual_sign -1 ddeltaU_ddisp sign -1;
41  jac += _test[_i][_qp] * jacsd * _vars[component_j]->phiFace()[_j][_qp];
42  jac -= _test[_i][_qp] * JacLD(component_j, /*neighbor=*/false);
43  break;
44  case Moose::ElementNeighbor: // Residual_sign -1 ddeltaU_ddisp sign 1;
45  jac -= _test[_i][_qp] * jacsd * _vars[component_j]->phiFaceNeighbor()[_j][_qp];
46  jac -= _test[_i][_qp] * JacLD(component_j, /*neighbor=*/true);
47  break;
48  case Moose::NeighborElement: // Residual_sign 1 ddeltaU_ddisp sign -1;
49  jac -= _test_neighbor[_i][_qp] * jacsd * _vars[component_j]->phiFace()[_j][_qp];
50  jac += _test_neighbor[_i][_qp] * JacLD(component_j, /*neighbor=*/false);
51  break;
52  case Moose::NeighborNeighbor: // Residual_sign 1 ddeltaU_ddisp sign 1;
53  jac += _test_neighbor[_i][_qp] * jacsd * _vars[component_j]->phiFaceNeighbor()[_j][_qp];
54  jac += _test_neighbor[_i][_qp] * JacLD(component_j, /*neighbor=*/true);
55  break;
56  }
57 
58  return jac;
59 }
60 
61 Real
62 CZMInterfaceKernelTotalLagrangian::JacLD(const unsigned int cc, const bool neighbor) const
63 {
64  Real jacld = 0;
65  RealVectorValue phi;
66  if (neighbor)
67  phi = 0.5 * _vars[cc]->gradPhiFaceNeighbor()[_j][_qp];
68  else
69  phi = 0.5 * _vars[cc]->gradPhiFace()[_j][_qp];
70 
71  for (unsigned int j = 0; j < 3; j++)
72  jacld += _dPK1traction_dF[_qp](_component, cc, j) * phi(j);
73 
74  mooseAssert(std::isfinite(jacld), "CZMInterfaceKernelTotalLagrangian JacLD is not finite");
75  return jacld;
76 }
NeighborElement
registerMooseObject("SolidMechanicsApp", CZMInterfaceKernelTotalLagrangian)
T & set(const std::string &name, bool quiet_mode=false)
Base class for implementing DG cohesive zone models (CZM) for 1D,2D, and 3D traction separation laws...
static InputParameters validParams()
CZMInterfaceKernelTotalLagrangian(const InputParameters &parameters)
ElementElement
Real computeDResidualDDisplacement(const unsigned int &component_j, const Moose::DGJacobianType &type) const override
method computing the derivative of residual[_component] w.r.t displacement[component_j] ...
DG cohesive zone model kernel for the Total Lagrangian formulation.
const unsigned int _component
the displacement component this kernel is operating on (0=x, 1=y, 2 =z)
std::vector< MooseVariable * > _vars
ElementNeighbor
const MaterialProperty< RankThreeTensor > & _dPK1traction_dF
Real JacLD(const unsigned int cc, const bool neighbor) const
method computing the jacobian contribution due to rotations and area changes
DGJacobianType
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const MaterialProperty< RankTwoTensor > & _dtraction_djump_global
NeighborNeighbor