https://mooseframework.inl.gov
ADArrayNodalKernel.h
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 #pragma once
11 
12 #include "NodalKernelBase.h"
13 #include "MooseVariableInterface.h"
14 
18 class ADArrayNodalKernel : public NodalKernelBase, public MooseVariableInterface<RealEigenVector>
19 {
20 public:
26 
28 
35  virtual void computeResidual() override;
36 
43  virtual void computeJacobian() override;
44 
51  virtual void computeOffDiagJacobian(unsigned int jvar) override;
52 
57  const MooseVariableFE<RealEigenVector> & variable() const override { return _var; }
58 
59  virtual void jacobianSetup() override { _my_node = nullptr; }
60 
61 protected:
65  virtual void computeQpResidual(ADRealEigenVector & residual) = 0;
66 
70  {
71  mooseError("I'm an AD object, so computeQpJacobian should never be called");
72  }
73 
76 
79 
81  const unsigned int _count;
82 
83 private:
86 
88  const Node * _my_node = nullptr;
89 };
Eigen::Matrix< ADReal, Eigen::Dynamic, 1 > ADRealEigenVector
Definition: MooseTypes.h:148
virtual void computeOffDiagJacobian(unsigned int jvar) override
Compute the off-diagonal Jacobian at one node.
ADRealEigenVector _work_vector
Work vector for residual.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
virtual void computeJacobian() override
Compute the Jacobian at one node.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
static InputParameters validParams()
Class constructor.
const MooseVariableFE< RealEigenVector > & variable() const override
Gets the variable this is active on.
virtual void computeResidual() override
Compute the residual at the current node.
virtual RealEigenVector computeQpJacobian()
Dummy method so we can make derived generic classes that template on <bool is_ad>=""> ...
const unsigned int _count
Number of components of the array variable.
Base class for creating nodal kernels with hand-coded Jacobians.
const ADArrayVariableValue & _u
Value of the unknown variable this is acting on.
MooseVariableFE< RealEigenVector > & _var
variable this works on
const Node * _my_node
Cache variable to make sure we don&#39;t do duplicate AD computations.
virtual void jacobianSetup() override
Gets called just before the Jacobian is computed and before this object is asked to do its job...
forward declarations
Base class for creating new types of nodal kernels.
Interface for objects that need to get values of MooseVariables.
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:271
ADArrayNodalKernel(const InputParameters &parameters)
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
virtual void computeQpResidual(ADRealEigenVector &residual)=0
The user can override this function to compute the residual at a node.