Line data Source code
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 : 16 : InputParameters 17 117967 : FVFluxBC::validParams() 18 : { 19 117967 : InputParameters params = FVBoundaryCondition::validParams(); 20 117967 : params += TwoMaterialPropertyInterface::validParams(); 21 117967 : params.registerSystemAttributeName("FVFluxBC"); 22 : 23 : // FVFluxBCs always rely on Boundary MaterialData 24 117967 : params.set<Moose::MaterialDataType>("_material_data_type") = Moose::BOUNDARY_MATERIAL_DATA; 25 117967 : params.set<bool>("_residual_object") = true; 26 : 27 117967 : return params; 28 0 : } 29 : 30 1955 : FVFluxBC::FVFluxBC(const InputParameters & parameters) 31 : : FVBoundaryCondition(parameters), 32 : NeighborCoupleableMooseVariableDependencyIntermediateInterface( 33 : this, /*nodal=*/false, /*neighbor_nodal=*/false, /*is_fv=*/true), 34 : TwoMaterialPropertyInterface(this, Moose::EMPTY_BLOCK_IDS, boundaryIDs()), 35 3910 : _u(_var.adSln()), 36 1955 : _u_neighbor(_var.adSlnNeighbor()) 37 : { 38 1955 : if (_var.kind() == Moose::VarKindType::VAR_AUXILIARY) 39 0 : paramError("variable", 40 : "There should not be a need to specify a flux " 41 : "boundary condition for an auxiliary variable."); 42 1955 : } 43 : 44 : void 45 481463 : FVFluxBC::updateCurrentFace(const FaceInfo & fi) 46 : { 47 481463 : _face_info = &fi; 48 481463 : _normal = fi.normal(); 49 481463 : _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. 58 481463 : if (_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) 59 0 : _normal = -_normal; 60 481463 : } 61 : 62 : void 63 481439 : FVFluxBC::computeResidual(const FaceInfo & fi) 64 : { 65 481439 : updateCurrentFace(fi); 66 : 67 481439 : if (_face_type == FaceInfo::VarFaceNeighbors::BOTH) 68 5 : mooseError("A FVFluxBC is being triggered on an internal face with centroid: ", 69 : fi.faceCentroid()); 70 481434 : else if (_face_type == FaceInfo::VarFaceNeighbors::NEITHER) 71 5 : 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 : 75 481429 : auto r = MetaPhysicL::raw_value(fi.faceArea() * fi.faceCoord() * computeQpResidual()); 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. 81 481429 : if (_face_type == FaceInfo::VarFaceNeighbors::ELEM) 82 481429 : prepareVectorTag(_assembly, _var.number()); 83 0 : else if (_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) 84 0 : prepareVectorTagNeighbor(_assembly, _var.number()); 85 : 86 481429 : _local_re(0) = r; 87 481429 : accumulateTaggedLocalResidual(); 88 481429 : } 89 : 90 : void 91 428 : FVFluxBC::computeResidualAndJacobian(const FaceInfo & fi) 92 : { 93 428 : computeJacobian(fi); 94 428 : } 95 : 96 : void 97 210268 : FVFluxBC::computeJacobian(const FaceInfo & fi) 98 : { 99 210268 : _face_info = &fi; 100 210268 : _normal = fi.normal(); 101 210268 : _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. 110 210268 : if (_face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR) 111 0 : _normal = -_normal; 112 : 113 210268 : ADReal r = fi.faceArea() * fi.faceCoord() * computeQpResidual(); 114 : 115 210268 : const auto & dof_indices = (_face_type == FaceInfo::VarFaceNeighbors::ELEM) 116 210268 : ? _var.dofIndices() 117 0 : : _var.dofIndicesNeighbor(); 118 : 119 : mooseAssert(dof_indices.size() == 1, "We're currently built to use CONSTANT MONOMIALS"); 120 : 121 420536 : addResidualsAndJacobian(_assembly, std::array<ADReal, 1>{{r}}, dof_indices, _var.scalingFactor()); 122 420536 : } 123 : 124 : const ADReal & 125 0 : FVFluxBC::uOnUSub() const 126 : { 127 : mooseAssert(_face_info, "The face info has not been set"); 128 0 : const auto ft = _face_info->faceType(std::make_pair(_var.number(), _var.sys().number())); 129 : mooseAssert( 130 : ft == FaceInfo::VarFaceNeighbors::ELEM || ft == FaceInfo::VarFaceNeighbors::NEIGHBOR, 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 : 139 0 : if (ft == FaceInfo::VarFaceNeighbors::ELEM) 140 0 : return _u[_qp]; 141 : else 142 0 : return _u_neighbor[_qp]; 143 : } 144 : 145 : const ADReal & 146 0 : FVFluxBC::uOnGhost() const 147 : { 148 : mooseAssert(_face_info, "The face info has not been set"); 149 0 : const auto ft = _face_info->faceType(std::make_pair(_var.number(), _var.sys().number())); 150 : mooseAssert( 151 : ft == FaceInfo::VarFaceNeighbors::ELEM || ft == FaceInfo::VarFaceNeighbors::NEIGHBOR, 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 : 160 0 : if (ft == FaceInfo::VarFaceNeighbors::ELEM) 161 0 : return _u_neighbor[_qp]; 162 : else 163 0 : return _u[_qp]; 164 : } 165 : 166 : Moose::ElemArg 167 0 : FVFluxBC::elemArg(const bool correct_skewness) const 168 : { 169 0 : return {&_face_info->elem(), correct_skewness}; 170 : } 171 : 172 : Moose::ElemArg 173 0 : FVFluxBC::neighborArg(const bool correct_skewness) const 174 : { 175 0 : return {_face_info->neighborPtr(), correct_skewness}; 176 : }