www.mooseframework.org
ADNodalBC.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "MooseVariableInterface.h"
13 #include "NodalBCBase.h"
14 #include "ADDirichletBCBase.h"
15 #include "ADFunctorInterface.h"
16 
20 template <typename T, typename Base>
21 class ADNodalBCTempl : public Base, public MooseVariableInterface<T>, public ADFunctorInterface
22 {
23 public:
25 
26  ADNodalBCTempl(const InputParameters & parameters);
27 
28  const MooseVariableFE<T> & variable() const override { return _var; }
29 
30  bool shouldSetComp(unsigned short i) const { return _set_components[i]; }
31 
32 protected:
36  virtual typename Moose::ADType<T>::type computeQpResidual() = 0;
37 
40 
42  const Node * const & _current_node;
43 
45  const unsigned int _qp = 0;
46 
48  const typename Moose::ADType<T>::type & _u;
49 
50  const std::array<bool, 3> _set_components;
51 
52  using Base::_fe_problem;
53  using Base::_subproblem;
54  using Base::_sys;
55  using Base::_tid;
56  using Base::addJacobian;
57  using Base::addMooseVariableDependency;
58  using Base::setResidual;
59 
60 private:
61  void computeResidual() override final;
62  void computeJacobian() override final;
63  void computeResidualAndJacobian() override;
64  void computeOffDiagJacobian(unsigned int jvar) override final;
65  void computeOffDiagJacobianScalar(unsigned int jvar) override final;
66 
70  template <typename ADResidual>
71  void addResidual(const ADResidual & residual, const std::vector<dof_id_type> & dof_indices);
72 
76  template <typename ADResidual>
77  void addJacobian(const ADResidual & residual, const std::vector<dof_id_type> & dof_indices);
78 
82 };
83 
84 template <>
86 
87 template <>
89 
void addJacobian(const ADResidual &residual, const std::vector< dof_id_type > &dof_indices)
process the Jacobian into the global data structures
Definition: ADNodalBC.C:113
Keeps track of stuff related to assembling.
Definition: Assembly.h:93
void computeResidualAndJacobian() override
Definition: ADNodalBC.C:157
Base class for deriving any automatic differentiation boundary condition of a integrated type...
Definition: ADNodalBC.h:21
Base class for automatic differentiation Dirichlet BCs.
void computeJacobian() override final
Definition: ADNodalBC.C:144
ADNodalBCTempl(const InputParameters &parameters)
Definition: ADNodalBC.C:55
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const Moose::ADType< T >::type & _u
Value of the unknown variable this BC is acting on.
Definition: ADNodalBC.h:48
const std::array< bool, 3 > _set_components
Definition: ADNodalBC.h:50
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
virtual Moose::ADType< T >::type computeQpResidual()=0
Compute this NodalBC&#39;s contribution to the residual at the current quadrature point.
Assembly & _undisplaced_assembly
A reference to the undisplaced assembly in order to ensure data gets correctly incorporated into the ...
Definition: ADNodalBC.h:81
void computeOffDiagJacobian(unsigned int jvar) override final
Definition: ADNodalBC.C:171
static InputParameters validParams()
Definition: ADNodalBC.C:25
void computeResidual() override final
Definition: ADNodalBC.C:131
MooseVariableFE< T > & _var
The variable that this NodalBC operates on.
Definition: ADNodalBC.h:39
Base class for deriving any boundary condition that works at nodes.
Definition: NodalBCBase.h:26
const unsigned int _qp
Pseudo-"quadrature point" index (Always zero for the current node)
Definition: ADNodalBC.h:45
const Node *const & _current_node
current node being processed
Definition: ADNodalBC.h:42
bool shouldSetComp(unsigned short i) const
Definition: ADNodalBC.h:30
void computeOffDiagJacobianScalar(unsigned int jvar) override final
Definition: ADNodalBC.C:180
Interface for objects that need to get values of MooseVariables.
const MooseVariableFE< T > & variable() const override
Definition: ADNodalBC.h:28
void addResidual(const ADResidual &residual, const std::vector< dof_id_type > &dof_indices)
process the residual into the global data structures
Definition: ADNodalBC.C:99
uint8_t dof_id_type