https://mooseframework.inl.gov
WeightedGapAux.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 
10 #include "WeightedGapAux.h"
11 
13 
16 {
18  params.addClassDescription(
19  "Returns the specified variable as an auxiliary variable with the same value.");
20  params.set<bool>("interpolate_normals") = false;
21  return params;
22 }
23 
25  : MortarNodalAuxKernel(parameters), _weighted_gap(0), _qp_gap_nodal(0)
26 {
27 }
28 
29 Real
31 {
32  _weighted_gap = 0;
33  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
34  {
36  for (_i = 0; _i < _test_lower.size(); ++_i)
38  }
39 
40  return _weighted_gap;
41 }
42 
43 void
45 {
46  const auto gap_vec = _phys_points_primary[_qp] - _phys_points_secondary[_qp];
47 
48  _qp_gap_nodal = gap_vec * (_JxW_msm[_qp] * _coord_msm[_qp]);
49 
51 }
52 
53 void
55 {
57 }
Real _msm_volume
The mortar segment volume.
static InputParameters validParams()
WeightedGapAux(const InputParameters &parameters)
Factory constructor, takes parameters so that all derived classes can be built using the same constru...
void computeQpProperties()
static InputParameters validParams()
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const MooseArray< Point > & _phys_points_primary
The locations of the quadrature points on the interior primary elements.
const std::vector< Real > & _JxW_msm
The element Jacobian times weights.
const libMesh::QBase *const & _qrule_msm
The quadrature rule on the mortar segment element.
Compute a weighted gap based on a mortar discretization.
const OutputTools< ComputeValueType >::VariableTestValue & _test_lower
The shape functions for the variable associated with the lower-dimensional element, e.g.
unsigned int _i
The current test function index.
Real _weighted_gap
The weighted gap.
std::vector< Point > _normals
the normals
const MooseArray< Real > & _coord_msm
Member for handling change of coordinate systems (xyz, rz, spherical) on mortar elements.
unsigned int n_points() const
const MooseArray< Point > & _phys_points_secondary
The locations of the quadrature points on the interior secondary elements.
Base class for creating new nodally-based mortar auxiliary kernels.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
RealVectorValue _qp_gap_nodal
The gap vector at the current quadrature point, used when we are not interpolating the normal vector...
void computeQpIProperties()
Real computeValue() override
Compute and return the value of the aux variable.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
unsigned int _qp
The current quadrature point index.
registerMooseObject("MooseApp", WeightedGapAux)