https://mooseframework.inl.gov
ADNodalKernel.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 "ADNodalKernel.h"
11 #include "SubProblem.h"
12 #include "SystemBase.h"
13 #include "MooseVariableFE.h"
14 #include "Assembly.h"
15 
16 #include "metaphysicl/raw_type.h"
17 
20 {
21  auto params = NodalKernelBase::validParams();
23  return params;
24 }
25 
27  : NodalKernelBase(parameters),
28  ADFunctorInterface(this),
30  true,
31  "variable",
34  _var(*mooseVariable()),
35  _u(_var.adDofValues())
36 {
38 }
39 
40 void
42 {
43  if (_var.isNodalDefined())
44  {
45  const auto dof_idx = _var.nodalDofIndex();
46  _qp = 0;
49  std::array<Real, 1>{{res}},
50  std::array<dof_id_type, 1>{{dof_idx}},
52  }
53 }
54 
55 void
57 {
58  if (_var.isNodalDefined())
59  {
60  const auto dof_idx = _var.nodalDofIndex();
61  _qp = 0;
62  const auto res = computeQpResidual();
64  std::array<ADReal, 1>{{res}},
65  std::array<dof_id_type, 1>{{dof_idx}},
67  }
68 }
69 
70 void
71 ADNodalKernel::computeOffDiagJacobian(const unsigned int jvar)
72 {
73  if (jvar == _var.number())
75 }
VarFieldType
Definition: MooseTypes.h:770
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.
static InputParameters validParams()
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:100
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
MooseVariable & _var
variable this works on
Definition: ADNodalKernel.h:67
MooseVariableFE< Real > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
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.
virtual ADReal computeQpResidual()=0
The user can override this function to compute the residual at a node.
void computeOffDiagJacobian(unsigned int jvar) override final
This method simply routes to computeJacobian whenever jvar == _var.number() since global AD computes ...
Definition: ADNodalKernel.C:71
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
static InputParameters validParams()
Class constructor.
Definition: ADNodalKernel.C:19
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
static InputParameters validParams()
Class constructor.
Base class for creating new types of nodal kernels.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual bool isNodalDefined() const override
Is this variable defined at nodes.
const dof_id_type & nodalDofIndex() const override
Interface for objects that need to get values of MooseVariables.
ADNodalKernel(const InputParameters &parameters)
Definition: ADNodalKernel.C:26
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
void computeJacobian() override
Compute the Jacobian at one node.
Definition: ADNodalKernel.C:56
unsigned int _qp
Quadrature point index.
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
void computeResidual() override
Compute the residual at the current node.
Definition: ADNodalKernel.C:41