https://mooseframework.inl.gov
ADIntegratedBC.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 "IntegratedBCBase.h"
13 #include "MooseVariableInterface.h"
14 #include "ADFunctorInterface.h"
15 
19 template <typename T>
21  public MooseVariableInterface<T>,
22  public ADFunctorInterface
23 {
24 public:
26 
28 
29  const MooseVariableFE<T> & variable() const override { return _var; }
30 
31 protected:
32  void computeResidual() override;
33  void computeJacobian() override;
34  void computeResidualAndJacobian() override;
35  void computeOffDiagJacobian(unsigned int jvar) override;
36  void computeOffDiagJacobianScalar(unsigned int jvar) override;
37 
45  virtual void computeResidualsForJacobian();
46 
50  virtual ADReal computeQpResidual() = 0;
51 
54 
57 
60 
61  // test functions
62 
67 
72 
75 
78 
81 
83  const bool _use_displaced_mesh;
84 
87  std::vector<Real> _residuals;
88  std::vector<ADReal> _residuals_and_jacobians;
89 
90 private:
94  void computeADJacobian();
95 };
96 
const ADTemplateVariableValue< T > & _u
the values of the unknown variable this BC is acting on
const bool _use_displaced_mesh
Whether this object is acting on the displaced mesh.
typename OutputTools< T >::VariableTestValue ADTemplateVariableTestValue
Definition: MooseTypes.h:623
void computeResidual() override
Compute this object&#39;s contribution to the residual.
typename OutputTools< T >::VariablePhiValue ADTemplateVariablePhiValue
Definition: MooseTypes.h:625
typename OutputTools< typename Moose::ADType< T >::type >::VariableValue ADTemplateVariableValue
Definition: MooseTypes.h:603
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
std::vector< ADReal > _residuals_and_jacobians
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:47
static InputParameters validParams()
const MooseArray< ADPoint > & _ad_q_points
(physical) quadrature points
const ADTemplateVariableGradient< T > & _grad_u
the gradient of the unknown variable this BC is acting on
An interface for accessing Moose::Functors for systems that care about automatic differentiation, e.g.
ADIntegratedBCTempl(const InputParameters &parameters)
void computeADJacobian()
compute all the Jacobian entries
typename OutputTools< typename Moose::ADType< T >::type >::VariableTestGradient ADTemplateVariableTestGradient
Definition: MooseTypes.h:631
const MooseArray< ADReal > & _ad_coord
The AD version of coord.
typename OutputTools< typename Moose::ADType< T >::type >::VariableGradient ADTemplateVariableGradient
Definition: MooseTypes.h:606
Base class for deriving any boundary condition of a integrated type.
MooseVariableFE< T > & _var
The variable that this IntegratedBC operates on.
const ADTemplateVariablePhiValue< T > & _phi
The current shape functions.
const ADTemplateVariableTestGradient< T > & _grad_test
gradients of test functions (in QPs)
void computeOffDiagJacobian(unsigned int jvar) override
Computes this object&#39;s contribution to off-diagonal blocks of the system Jacobian matrix...
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
forward declarations
const ADTemplateVariableTestValue< T > & _test
test function values (in QPs)
ADReal _r
Data members for holding residuals.
Interface for objects that need to get values of MooseVariables.
std::vector< Real > _residuals
Base class for deriving any boundary condition of a integrated type.
const InputParameters & parameters() const
Get the parameters of the object.
const MooseArray< ADReal > & _ad_JxW
The ad version of JxW.
const MooseVariableFE< T > & variable() const override
Returns the variable that this object operates on.
virtual void computeResidualsForJacobian()
compute the _residuals member for filling the Jacobian.
const MooseArray< ADPoint > & _normals
normals at quadrature points
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
virtual ADReal computeQpResidual()=0
Compute this IntegratedBC&#39;s contribution to the residual at the current quadrature point...
void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.