https://mooseframework.inl.gov
ScalarLMKernel.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 "ScalarLMKernel.h"
11 
14 
15 template <bool is_ad>
18 {
20  params.addClassDescription("This class is used to enforce integral of phi = V_0 with a "
21  "Lagrange multiplier approach.");
22  params.renameCoupledVar("scalar_variable", "kappa", "Primary coupled scalar variable");
23  params.addRequiredParam<PostprocessorName>(
24  "pp_name", "Name of the Postprocessor containing the volume of the domain.");
25  params.addRequiredParam<Real>(
26  "value", "Given (constant) which we want the integral of the solution variable to match.");
27 
28  return params;
29 }
30 
31 template <bool is_ad>
33  : GenericKernelScalar<is_ad>(parameters),
34  _value(this->template getParam<Real>("value")),
35  _pp_value(this->getPostprocessorValue("pp_name"))
36 {
37 }
38 
39 template <bool is_ad>
42 {
43  return _kappa[0] * _test[_i][_qp];
44 }
45 
46 template <bool is_ad>
49 {
50  return _u[_qp] - _value / _pp_value;
51 }
52 
53 template <bool is_ad>
54 Real
56 {
57  // This function will never be called for the AD version. But because C++ does
58  // not support an optional function declaration based on a template parameter,
59  // we must keep this template for all cases.
60  mooseAssert(
61  !is_ad,
62  "In ADScalarLMKernel, computeScalarQpJacobian should not be called. Check computeJacobian "
63  "implementation.");
64  return 0.;
65 }
66 
67 template <bool is_ad>
68 Real
70 {
71  // This function will never be called for the AD version. But because C++ does
72  // not support an optional function declaration based on a template parameter,
73  // we must keep this template for all cases.
74  mooseAssert(!is_ad,
75  "In ADScalarLMKernel, computeQpOffDiagJacobianScalar should not be called. Check "
76  "computeOffDiagJacobianScalar "
77  "implementation.");
78  if (svar == _kappa_var)
79  return _test[_i][_qp];
80  else
81  return 0.;
82 }
83 
84 template <bool is_ad>
85 Real
87 {
88  // This function will never be called for the AD version. But because C++ does
89  // not support an optional function declaration based on a template parameter,
90  // we must keep this template for all cases.
91  mooseAssert(!is_ad,
92  "In ADScalarLMKernel, computeScalarQpOffDiagJacobian should not be called. Check "
93  "computeOffDiagJacobian "
94  "implementation.");
95  if (jvar == _var.number())
96  return _phi[_j][_qp];
97  else
98  return 0.;
99 }
100 
101 template class ScalarLMKernelTempl<false>;
102 template class ScalarLMKernelTempl<true>;
Moose::GenericType< Real, is_ad > GenericReal
Definition: MooseTypes.h:648
virtual GenericReal< is_ad > computeQpResidual() override
Method for computing the residual at quadrature points.
virtual Real computeScalarQpOffDiagJacobian(const unsigned int jvar) override
Method for computing d-_kappa-residual / d-_var at quadrature points.
static InputParameters validParams()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
void renameCoupledVar(const std::string &old_name, const std::string &new_name, const std::string &new_docstring)
Rename a coupled variable and provide a new documentation string.
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...
virtual Real computeQpOffDiagJacobianScalar(const unsigned int svar) override
Method for computing d-_var-residual / d-_kappa at quadrature points.
virtual GenericReal< is_ad > computeScalarQpResidual() override
Method for computing the scalar part of residual at quadrature points.
This Kernel implements part of the equation that enforces the constraint of.
virtual Real computeScalarQpJacobian() override
Method for computing the scalar variable part of Jacobian at quadrature points.
registerMooseObject("MooseApp", ScalarLMKernel)
ScalarLMKernelTempl(const InputParameters &parameters)
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()