https://mooseframework.inl.gov
FVScalarLagrangeMultiplierInterface.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 #include "MooseVariableScalar.h"
12 #include "Assembly.h"
13 
16 {
18  params.addClassDescription(
19  "This class should be inherited to create interface penalties in finite volume.");
20  params.addRequiredCoupledVar("lambda", "The name of the scalar lagrange multiplier");
21 
22  return params;
23 }
24 
26  const InputParameters & params)
27  : FVInterfaceKernel(params),
28  _lambda_var(*getScalarVar("lambda", 0)),
29  _lambda(adCoupledScalarValue("lambda"))
30 {
31  if (var1().sys().number() != var2().sys().number())
32  mooseError(this->type(), " does not support multiple nonlinear systems!");
33 }
34 
35 void
37 {
38  setupData(fi);
39 
40  const auto var_elem_num = _elem_is_one ? var1().number() : var2().number();
41  const auto var_neigh_num = _elem_is_one ? var2().number() : var1().number();
42 
43  const auto r =
44  MetaPhysicL::raw_value(_lambda[0]) * fi.faceArea() * fi.faceCoord() * (2 * _elem_is_one - 1);
45 
46  // Primal residual
47  addResidual(r, var_elem_num, false);
48  addResidual(-r, var_neigh_num, true);
49 
50  // LM residual. We may not have any actual ScalarKernels in our simulation so we need to manually
51  // make sure the scalar residuals get cached for later addition
52  const auto lm_r = MetaPhysicL::raw_value(computeQpResidual()) * fi.faceArea() * fi.faceCoord();
54  std::array<Real, 1>{{lm_r}},
57 }
58 
59 void
61 {
62  setupData(fi);
63 
64  const auto & elem_dof_indices = _elem_is_one ? var1().dofIndices() : var2().dofIndices();
65  const auto & neigh_dof_indices =
67  mooseAssert((elem_dof_indices.size() == 1) && (neigh_dof_indices.size() == 1),
68  "We're currently built to use CONSTANT MONOMIALS");
69  const auto elem_scaling_factor = _elem_is_one ? var1().scalingFactor() : var2().scalingFactor();
70  const auto neigh_scaling_factor = _elem_is_one ? var2().scalingFactor() : var1().scalingFactor();
71 
72  // Primal
73  const auto primal_r = _lambda[0] * fi.faceArea() * fi.faceCoord() * (2 * _elem_is_one - 1);
75  _assembly, std::array<ADReal, 1>{{primal_r}}, elem_dof_indices, elem_scaling_factor);
77  _assembly, std::array<ADReal, 1>{{-primal_r}}, neigh_dof_indices, neigh_scaling_factor);
78 
79  // LM
80  const auto lm_r = computeQpResidual() * fi.faceArea() * fi.faceCoord();
81  mooseAssert(_lambda_var.dofIndices().size() == 1, "We should only have one dof");
83  std::array<ADReal, 1>{{lm_r}},
86 }
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
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...
const MooseVariableFV< Real > & var1() const
unsigned int number() const
Get variable number coming from libMesh.
void computeResidual(const FaceInfo &fi) override final
Compute the residual on the supplied face.
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...
Real & faceCoord()
Sets/gets the coordinate transformation factor (for e.g.
Definition: FaceInfo.h:64
static InputParameters validParams()
Real faceArea() const
Returns the face area of face id.
Definition: FaceInfo.h:60
bool _elem_is_one
Whether the current element is associated with variable/subdomain 1 or 2.
Assembly & _assembly
The Assembly object.
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
Base class for creating kernels that interface physics between subdomains.
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
const ADVariableValue & _lambda
The Lagrange Multiplier value.
void setupData(const FaceInfo &fi)
setup data useful for this object
void computeJacobian(const FaceInfo &fi) override final
Compute the jacobian on the supplied face.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
const MooseVariableScalar & _lambda_var
The Lagrange Multiplier variable.
FVScalarLagrangeMultiplierInterface(const InputParameters &params)
ADReal computeQpResidual() override=0
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:267
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 addResidual(Real resid, unsigned int var_num, bool neighbor)
Process the provided residual given var_num and whether this is on the neighbor side.
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
const MooseVariableFV< Real > & var2() const
virtual const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.