https://mooseframework.inl.gov
ADDirichletBCBaseTempl.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 "ADDirichletBCBaseTempl.h"
11 #include "MooseVariableFE.h"
12 #include "SystemBase.h"
13 #include "libmesh/node.h"
14 
15 template <typename T>
18 {
20  params.addParam<bool>(
21  "preset", true, "Whether or not to preset the BC (apply the value before the solve begins).");
22  return params;
23 }
24 
25 template <typename T>
27  : ADNodalBCTempl<T, ADDirichletBCBase>(parameters)
28 {
29 }
30 
31 template <typename T>
32 void
34 {
35  mooseAssert(this->_preset, "BC is not preset");
36 
37  if (_var.isNodalDefined())
38  {
39  const auto n_comp = _current_node->n_comp(_sys.number(), _var.number());
40  const auto value = MetaPhysicL::raw_value(computeQpValue());
41  for (const auto i : make_range(n_comp))
42  {
43  const auto dof_idx = _current_node->dof_number(_sys.number(), _var.number(), i);
44  if constexpr (std::is_same<T, Real>::value)
45  {
46  mooseAssert(n_comp == 1, "This should only be unity");
47  current_solution.set(dof_idx, value);
48  }
49  else
50  {
51  if (shouldSetComp(i))
52  current_solution.set(dof_idx, value(i));
53  }
54  }
55  }
56 }
57 
58 template <typename T>
61 {
62  return _u - computeQpValue();
63 }
64 
65 template class ADDirichletBCBaseTempl<Real>;
virtual void computeValue(NumericVector< Number > &current_solution) override
Method to preset the nodal value if applicable.
Base class for deriving any automatic differentiation boundary condition of a integrated type...
Definition: ADNodalBC.h:21
Base class for automatic differentiation Dirichlet BCs.
virtual Moose::ADType< T >::type computeQpResidual() override
Compute this NodalBC&#39;s contribution to the residual at the current quadrature point.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
ADDirichletBCBaseTempl(const InputParameters &parameters)
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
static InputParameters validParams()
Definition: ADNodalBC.C:26
IntRange< T > make_range(T beg, T end)
static InputParameters validParams()
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
virtual void set(const numeric_index_type i, const Number value)=0