https://mooseframework.inl.gov
OneDEqualValueConstraintBC.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("Computes the integral of lambda times dg term from the mortar method"
19  " (for two 1D domains only).");
20  params.addRequiredCoupledVar("lambda", "Lagrange multiplier");
21  params.addRequiredParam<unsigned int>("component", "Component of the Lagrange multiplier");
22  params.addRequiredParam<Real>(
23  "vg",
24  "Variation of the constraint g wrt this surface (+1 or -1). Note: g = value1 - value2 = 0 ");
25  return params;
26 }
27 
29  : IntegratedBC(parameters),
30  _lambda(coupledScalarValue("lambda")),
31  _lambda_var_number(coupledScalar("lambda")),
32  _component(getParam<unsigned int>("component")),
33  _vg(getParam<Real>("vg"))
34 {
35 }
36 
37 Real
39 {
40  return _lambda[_component] * _vg * _test[_i][_qp];
41 }
42 
43 Real
45 {
46  return 0.;
47 }
48 
49 Real
51 {
52  if (jvar == _lambda_var_number)
53  {
54  if (_j == _component)
55  return _vg * _test[_i][_qp];
56  else
57  return 0.;
58  }
59  else
60  return 0.;
61 }
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
virtual Real computeQpResidual() override
Method for computing the residual at quadrature points.
registerMooseObject("MooseApp", OneDEqualValueConstraintBC)
static InputParameters validParams()
Definition: IntegratedBC.C:23
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _i
i-th, j-th index for enumerating test and shape functions
This is the term from the mortar method.
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
unsigned int _qp
quadrature point index
virtual Real computeQpJacobian() override
Method for computing the diagonal Jacobian at quadrature points.
virtual Real computeQpOffDiagJacobianScalar(unsigned int jvar) override
Method for computing an off-diagonal jacobian component from a scalar var.
OneDEqualValueConstraintBC(const InputParameters &parameters)
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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...
static InputParameters validParams()
void ErrorVector unsigned int