https://mooseframework.inl.gov
ScalarLagrangeMultiplier.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 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseVariableFE.h"
15 #include "MooseVariableScalar.h"
16 #include "SystemBase.h"
17 
18 #include "libmesh/quadrature.h"
19 
21 
24 {
26  params.addClassDescription("This class is used to enforce integral of phi = V_0 with a "
27  "Lagrange multiplier approach.");
28  params.addRequiredCoupledVar("lambda", "Lagrange multiplier variable");
29 
30  return params;
31 }
32 
34  : Kernel(parameters), _lambda_var(coupledScalar("lambda")), _lambda(coupledScalarValue("lambda"))
35 {
36 }
37 
39 
40 Real
42 {
43  return _lambda[0] * _test[_i][_qp];
44 }
45 
46 void
48 {
49  // Note: Here we are assembling the contributions for _both_ Eq. (9) and (10)
50  // in the detailed writeup on this problem [0]. This is important to highlight
51  // because it is slightly different from the way things usually work in MOOSE
52  // because we are computing Jacobian contributions _for a different equation_
53  // than the equation which this Kernel is assigned in the input file. We do
54  // this for both simplicity (results in slightly less code) and efficiency:
55  // the contribution is symmetric, so it can be computed once and used twice.
56  //
57  // [0]: https://github.com/idaholab/large_media/blob/master/framework/scalar_constraint_kernel.pdf
58 
60 
62  for (_i = 0; _i < _test.size(); _i++)
63  for (_j = 0; _j < jv.order(); _j++)
64  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
67 
68  const auto ke_copy = _local_ke;
70  ke_copy.get_transpose(_local_ke);
72 }
73 
74 Real
76 {
77  if (jvar == _lambda_var)
78  return _test[_i][_qp];
79  else
80  return 0.;
81 }
static InputParameters validParams()
Definition: Kernel.C:24
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:72
virtual Real computeQpOffDiagJacobianScalar(unsigned int jvar) override
For coupling scalar variables.
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:52
unsigned int number() const
Get variable number coming from libMesh.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual Real computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
THREAD_ID _tid
The thread ID for this kernel.
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
SystemBase & _sys
Reference to the EquationSystem object.
const VariableValue & _lambda
Lagrange multiplier variable value.
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:49
unsigned int _i
current index for the test function
Definition: KernelBase.h:58
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
registerMooseObject("MooseApp", ScalarLagrangeMultiplier)
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz, spherical, etc.)
Definition: KernelBase.h:55
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
libMesh::Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
unsigned int _j
current index for the shape function
Definition: KernelBase.h:61
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:144
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Class for scalar variables (they are different).
Definition: Kernel.h:15
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...
ScalarLagrangeMultiplier(const InputParameters &parameters)
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the active tags.
unsigned int _lambda_var
Lagrange multiplier variable ID.
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43
This Kernel implements part of the equation that enforces the constraint of.