https://mooseframework.inl.gov
FVScalarLagrangeMultiplierConstraint.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 #include "MooseVariableScalar.h"
13 #include "MooseVariableFV.h"
14 #include "Assembly.h"
15 
18 {
20  params.addClassDescription(
21  "Base class for imposing constraints using scalar Lagrange multipliers");
22  params.addParam<PostprocessorName>("phi0", "0", "The value that the constraint will enforce.");
23  params.addRequiredCoupledVar("lambda", "Lagrange multiplier variable");
24  return params;
25 }
26 
28  const InputParameters & parameters)
29  : FVElementalKernel(parameters),
30  _phi0(getPostprocessorValue("phi0")),
31  _lambda_var(*getScalarVar("lambda", 0)),
32  _lambda(adCoupledScalarValue("lambda"))
33 {
34 }
35 
36 void
38 {
39  const auto volume = _assembly.elemVolume();
41  std::array<ADReal, 1>{{_lambda[0] * volume}},
42  _var.dofIndices(),
45  std::array<ADReal, 1>{{computeQpResidual() * volume}},
48 }
49 
50 void
52 {
53  // Primal residual
55  mooseAssert(_local_re.size() == 1, "We should only have a single dof");
56  mooseAssert(_lambda.size() == 1 && _lambda_var.order() == 1,
57  "The lambda variable should be first order");
60 
61  // LM residual. We may not have any actual ScalarKernels in our simulation so we need to manually
62  // make sure the scalar residuals get cached for later addition
64  mooseAssert(_lambda_var.dofIndices().size() == 1, "We should only have a single dof");
66  std::array<Real, 1>{{lm_r}},
69 }
70 
71 void
73 {
74 }
75 
76 void
78 {
79  // Primal
80  mooseAssert(_lambda.size() == 1 && _lambda_var.order() == 1,
81  "The lambda variable should be first order");
82  const auto primal_r = _lambda[0] * _assembly.elemVolume();
83  mooseAssert(_var.dofIndices().size() == 1, "We should only have one dof");
85  _assembly, std::array<ADReal, 1>{{primal_r}}, _var.dofIndices(), _var.scalingFactor());
86 
87  // LM
88  const auto lm_r = computeQpResidual() * _assembly.elemVolume();
89  mooseAssert(_lambda_var.dofIndices().size() == 1, "We should only have one dof");
91  std::array<ADReal, 1>{{lm_r}},
94 }
ADReal computeQpResidual() override=0
This is the primary function that must be implemented for flux kernel terms.
void addResidualsAndJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the provided d...
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
unsigned int number() const
Get variable number coming from libMesh.
void addResiduals(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals corresponding to the provided dof indices.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void computeJacobian() override final
Compute this object&#39;s contribution to the diagonal Jacobian entries.
const ADVariableValue & _lambda
The Lagrange Multiplier value.
void computeResidualAndJacobian() override final
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided residual derivatives into the Jacobian for the provided dof indices.
FVElemental is used for calculating residual contributions from volume integral terms of a PDE where ...
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
static InputParameters validParams()
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
void computeResidual() override final
Usually you should not override these functions - they have some tricky stuff in them that you don&#39;t ...
Real volume(const MeshBase &mesh, unsigned int dim=libMesh::invalid_uint)
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...
const MooseVariableScalar & _lambda_var
The Lagrange Multiplier variable.
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
virtual unsigned int size() const override final
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...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
const Real & elemVolume() const
Returns the reference to the current element volume.
Definition: Assembly.h:401
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
FVScalarLagrangeMultiplierConstraint(const InputParameters &parameters)
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
virtual const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
MooseVariableFV< Real > & _var