https://mooseframework.inl.gov
FVFluxBC.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 "FVFluxBC.h"
11 #include "MooseVariableFV.h"
12 #include "SystemBase.h"
13 #include "Assembly.h"
14 #include "ADUtils.h"
15 
18 {
21  params.registerSystemAttributeName("FVFluxBC");
22 
23  // FVFluxBCs always rely on Boundary MaterialData
24  params.set<Moose::MaterialDataType>("_material_data_type") = Moose::BOUNDARY_MATERIAL_DATA;
25  params.set<bool>("_residual_object") = true;
26 
27  return params;
28 }
29 
31  : FVBoundaryCondition(parameters),
33  this, /*nodal=*/false, /*neighbor_nodal=*/false, /*is_fv=*/true),
34  TwoMaterialPropertyInterface(this, Moose::EMPTY_BLOCK_IDS, boundaryIDs()),
35  _u(_var.adSln()),
36  _u_neighbor(_var.adSlnNeighbor())
37 {
39  paramError("variable",
40  "There should not be a need to specify a flux "
41  "boundary condition for an auxiliary variable.");
42 }
43 
44 void
46 {
47  _face_info = &fi;
48  _normal = fi.normal();
49  _face_type = fi.faceType(std::make_pair(_var.number(), _var.sys().number()));
50 
51  // For FV flux kernels, the normal is always oriented outward from the lower-id
52  // element's perspective. But for BCs, there is only a residual
53  // contribution to one element (one side of the face). Because of this, we
54  // make an exception and orient the normal to point outward from whichever
55  // side of the face the BC's variable is defined on; we flip it if this
56  // variable is defined on the neighbor side of the face (instead of elem) since
57  // the FaceInfo normal polarity is always oriented with respect to the lower-id element.
59  _normal = -_normal;
60 }
61 
62 void
64 {
66 
68  mooseError("A FVFluxBC is being triggered on an internal face with centroid: ",
69  fi.faceCentroid());
71  mooseError("A FVFluxBC is being triggered on a face which does not connect to a block ",
72  "with the relevant finite volume variable. Its centroid: ",
73  fi.faceCentroid());
74 
76 
77  // This could be an "internal" boundary - one created by variable block
78  // restriction where the var is only defined on one side of the face. We
79  // need to make sure that we add the residual contribution to the correct
80  // side - the one where the variable is defined.
85 
86  _local_re(0) = r;
88 }
89 
90 void
92 {
93  computeJacobian(fi);
94 }
95 
96 void
98 {
99  _face_info = &fi;
100  _normal = fi.normal();
101  _face_type = fi.faceType(std::make_pair(_var.number(), _var.sys().number()));
102 
103  // For FV flux kernels, the normal is always oriented outward from the lower-id
104  // element's perspective. But for BCs, there is only a Jacobian
105  // contribution to one element (one side of the face). Because of this, we
106  // make an exception and orient the normal to point outward from whichever
107  // side of the face the BC's variable is defined on; we flip it if this
108  // variable is defined on the neighbor side of the face (instead of elem) since
109  // the FaceInfo normal polarity is always oriented with respect to the lower-id element.
111  _normal = -_normal;
112 
113  ADReal r = fi.faceArea() * fi.faceCoord() * computeQpResidual();
114 
115  const auto & dof_indices = (_face_type == FaceInfo::VarFaceNeighbors::ELEM)
116  ? _var.dofIndices()
118 
119  mooseAssert(dof_indices.size() == 1, "We're currently built to use CONSTANT MONOMIALS");
120 
121  addResidualsAndJacobian(_assembly, std::array<ADReal, 1>{{r}}, dof_indices, _var.scalingFactor());
122 }
123 
124 const ADReal &
126 {
127  mooseAssert(_face_info, "The face info has not been set");
128  const auto ft = _face_info->faceType(std::make_pair(_var.number(), _var.sys().number()));
129  mooseAssert(
131  "The variable " << _var.name()
132  << " should be defined on exactly one adjacent subdomain for FVFluxBC "
133  << this->name());
134  mooseAssert(_qp == 0,
135  "At the time of writing, we only support one quadrature point, which should "
136  "correspond to the location of the cell centroid. If that changes, we should "
137  "probably change the body of FVFluxBC::uOnUSub");
138 
140  return _u[_qp];
141  else
142  return _u_neighbor[_qp];
143 }
144 
145 const ADReal &
147 {
148  mooseAssert(_face_info, "The face info has not been set");
149  const auto ft = _face_info->faceType(std::make_pair(_var.number(), _var.sys().number()));
150  mooseAssert(
152  "The variable " << _var.name()
153  << " should be defined on exactly one adjacent subdomain for FVFluxBC "
154  << this->name());
155  mooseAssert(_qp == 0,
156  "At the time of writing, we only support one quadrature point, which should "
157  "correspond to the location of the cell centroid. If that changes, we should "
158  "probably change the body of FVFluxBC::uOnGhost");
159 
161  return _u_neighbor[_qp];
162  else
163  return _u[_qp];
164 }
165 
167 FVFluxBC::elemArg(const bool correct_skewness) const
168 {
169  return {&_face_info->elem(), correct_skewness};
170 }
171 
173 FVFluxBC::neighborArg(const bool correct_skewness) const
174 {
175  return {_face_info->neighborPtr(), correct_skewness};
176 }
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const final
Get neighbor DOF indices for currently selected element.
const FaceInfo * _face_info
Holds information for the face we are currently examining.
void computeResidualAndJacobian(const FaceInfo &fi) override
Compute the residual and Jacobian on the supplied face.
Definition: FVFluxBC.C:91
FaceInfo::VarFaceNeighbors _face_type
The variable face type.
Definition: FVFluxBC.h:73
void addResidualsAndJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the provided d...
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
Intermediate base class that ties together all the interfaces for getting MooseVariables with the Moo...
Moose::ElemArg neighborArg(bool correct_skewness=false) const
Definition: FVFluxBC.C:173
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
static InputParameters validParams()
Definition: FVFluxBC.C:17
unsigned int number() const
Get variable number coming from libMesh.
const ADReal & uOnUSub() const
Definition: FVFluxBC.C:125
Base class for creating new types of boundary conditions.
const Elem & elem() const
Definition: FaceInfo.h:81
MaterialDataType
MaterialData types.
Definition: MooseTypes.h:692
const ADVariableValue & _u_neighbor
Definition: FVFluxBC.h:47
const Point & faceCentroid() const
Returns the coordinates of the face centroid.
Definition: FaceInfo.h:71
T & set(const std::string &name, bool quiet_mode=false)
Returns a writable reference to the named parameters.
Assembly & _assembly
Reference to assembly.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
void registerSystemAttributeName(const std::string &value)
This method is used to define the MOOSE system name that is used by the TheWarehouse object for stori...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
MooseVariableFV< Real > & _var
Real & faceCoord()
Sets/gets the coordinate transformation factor (for e.g.
Definition: FaceInfo.h:64
DualNumber< Real, DNDerivativeType, true > ADReal
Definition: ADRealForward.h:46
FVFluxBC(const InputParameters &parameters)
Definition: FVFluxBC.C:30
Real faceArea() const
Returns the face area of face id.
Definition: FaceInfo.h:60
const unsigned int _qp
Definition: FVFluxBC.h:45
const ADVariableValue & _u
Definition: FVFluxBC.h:46
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
const std::string & name() const
Get the name of the class.
Definition: MooseBase.h:99
const Elem * neighborPtr() const
Definition: FaceInfo.h:84
void updateCurrentFace(const FaceInfo &fi)
Update internal structures (normal, face type, etc) for the given face.
Definition: FVFluxBC.C:45
const ADReal & uOnGhost() const
Definition: FVFluxBC.C:146
ADRealVectorValue _normal
Definition: FVFluxBC.h:49
A structure that is used to evaluate Moose functors logically at an element/cell center.
const Point & normal() const
Returns the unit normal vector for the face oriented outward from the face&#39;s elem element...
Definition: FaceInfo.h:68
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1149
void prepareVectorTagNeighbor(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags for DG and interface kernels...
const std::set< SubdomainID > EMPTY_BLOCK_IDS
Definition: MooseTypes.h:684
Moose::VarKindType kind() const
Kind of the variable (Nonlinear, Auxiliary, ...)
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
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:267
This interface is designed for DGKernel, InternalSideUserObject, InterfaceUserObject, where material properties on a side of both its primary side (face) and its secondary side (neighbor) all required.
Moose::ElemArg elemArg(bool correct_skewness=false) const
Definition: FVFluxBC.C:167
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
static InputParameters validParams()
SystemBase & sys()
Get the system this variable is part of.
void computeJacobian(const FaceInfo &fi) override
Compute the jacobian on the supplied face.
Definition: FVFluxBC.C:97
virtual ADReal computeQpResidual()=0
void computeResidual(const FaceInfo &fi) override
Compute the residual on the supplied face.
Definition: FVFluxBC.C:63
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
virtual const std::vector< dof_id_type > & dofIndices() const final
Get local DoF indices.
VarFaceNeighbors faceType(const std::pair< unsigned int, unsigned int > &var_sys) const
Returns which side(s) the given variable-system number pair is defined on for this face...
Definition: FaceInfo.h:225