https://mooseframework.inl.gov
ADArrayNodalKernel.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 "ADArrayNodalKernel.h"
11 #include "SubProblem.h"
12 #include "SystemBase.h"
13 #include "MooseVariableFE.h"
14 #include "Assembly.h"
15 
18 {
20 }
21 
23  : NodalKernelBase(parameters),
25  true,
26  "variable",
29  _var(*mooseVariable()),
30  _u(_var.adDofValues()),
31  _count(_var.count()),
32  _work_vector(_count)
33 {
35 }
36 
37 void
39 {
40  if (!_var.isNodalDefined())
41  return;
42  _qp = 0;
46  _var.dofIndices(),
48 }
49 
50 void
52 {
53  if (!_var.isNodalDefined())
54  return;
55  _qp = 0;
58 }
59 
60 void
62 {
63  if (_my_node != _current_node)
64  {
67  }
68 }
VarFieldType
Definition: MooseTypes.h:770
const std::vector< Real > & arrayScalingFactor() const
virtual void computeOffDiagJacobian(unsigned int jvar) override
Compute the off-diagonal Jacobian at one node.
ADRealEigenVector _work_vector
Work vector for residual.
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.
virtual void computeJacobian() override
Compute the Jacobian at one node.
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...
static InputParameters validParams()
Class constructor.
virtual void computeResidual() override
Compute the residual at the current node.
MooseVariableFE< RealEigenVector > * 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.
const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
MooseVariableFE< RealEigenVector > & _var
variable this works on
const Node * _my_node
Cache variable to make sure we don&#39;t do duplicate AD computations.
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.
virtual bool isNodalDefined() const override
Is this variable defined at nodes.
Interface for objects that need to get values of MooseVariables.
ADArrayNodalKernel(const InputParameters &parameters)
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
virtual void computeQpResidual(ADRealEigenVector &residual)=0
The user can override this function to compute the residual at a node.
const Node *const & _current_node
current node being processed
unsigned int _qp
Quadrature point index.