23 #include "libmesh/elem.h" 24 #include "libmesh/vector_value.h" 25 #include "libmesh/tensor_value.h" 26 #include "libmesh/dense_vector.h" 27 #include "libmesh/dense_matrix.h" 28 #include "libmesh/libmesh_common.h" 37 template <
typename ActionFunctor>
44 const auto coord_type = subproblem.
getCoordSystem(elem.subdomain_id());
69 const bool two_term_expansion,
70 const bool correct_skewness,
71 const Elem * elem_to_extrapolate_from,
75 bool elem_to_extrapolate_from_is_fi_elem;
76 std::tie(elem_to_extrapolate_from, elem_to_extrapolate_from_is_fi_elem) =
77 [
this, &fi, elem_to_extrapolate_from]() -> std::pair<const Elem *, bool>
79 if (elem_to_extrapolate_from)
80 return {elem_to_extrapolate_from, elem_to_extrapolate_from == &fi.
elem()};
83 const auto [elem_guaranteed_to_have_dofs,
85 elem_guaranteed_to_have_dofs_is_fi_elem] =
88 return {elem_guaranteed_to_have_dofs, elem_guaranteed_to_have_dofs_is_fi_elem};
94 const Point vector_to_face = elem_to_extrapolate_from_is_fi_elem
102 fi, two_term_expansion, correct_skewness, elem_to_extrapolate_from, time);
104 return boundary_value;
110 const bool correct_skewness)
const 112 VectorValue<ADReal> unc_grad;
115 const auto & cell_gradient =
adGradSln(&fi.
elem(), time, correct_skewness);
116 const auto & normal = fi.
normal();
117 unc_grad = cell_gradient - (cell_gradient * normal) * normal;
125 const VectorValue<ADReal> &
128 bool correct_skewness)
const 147 std::vector<std::pair<const FaceInfo *, bool>> ebf_faces;
151 VectorValue<ADReal> &
grad = *value_pointer;
153 bool volume_set =
false;
189 std::vector<VectorValue<Real>> ebf_grad_coeffs;
192 std::vector<const ADReal *> ebf_b;
195 std::vector<VectorValue<Real>> grad_ebf_coeffs;
197 VectorValue<ADReal> grad_b = 0;
201 std::vector<TensorValue<Real>> fdf_grad_centroid_coeffs;
203 const unsigned int lm_dim = LIBMESH_DIM;
205 auto action_functor = [&volume_set,
216 &fdf_grad_centroid_coeffs,
219 this](
const Elem & functor_elem,
220 const Elem *
const neighbor,
222 const Point & surface_vector,
224 const bool elem_has_info)
226 mooseAssert(fi,
"We need a FaceInfo for this action_functor");
227 mooseAssert(elem == &functor_elem,
228 "Just a sanity check that the element being passed in is the one we passed out.");
235 ebf_faces.push_back(std::make_pair(fi, fdf_face));
238 ebf_grad_coeffs.push_back(-1. * (elem_has_info
239 ? (fi->faceCentroid() - fi->elemCentroid())
240 : (fi->faceCentroid() - fi->neighborCentroid())));
241 ebf_b.push_back(&elem_value);
244 grad_ebf_coeffs.push_back(-surface_vector);
250 fdf_grad_centroid_coeffs.emplace_back();
251 auto & current_coeffs = fdf_grad_centroid_coeffs.back();
252 const auto normal = fi->normal();
256 auto & current_coeff = current_coeffs(i,
j);
257 current_coeff = normal(i) * normal(
j);
268 grad_b += surface_vector * elem_value;
282 "We've run out of face types");
293 this->
_subproblem, functor_elem.subdomain_id(), fi->elemCentroid(), coord);
294 volume = fi->elemVolume() * coord;
299 this->
_subproblem, neighbor->subdomain_id(), fi->neighborCentroid(), coord);
300 volume = fi->neighborVolume() * coord;
309 mooseAssert(volume_set &&
volume > 0,
"We should have set the volume");
313 if (coord_system == Moose::CoordinateSystemType::COORD_RZ)
316 grad_b(r_coord) -= elem_value / elem->vertex_average()(r_coord);
320 coord_system != Moose::CoordinateSystemType::COORD_RSPHERICAL,
321 "We have not yet implemented the correct translation from gradient to divergence for " 322 "spherical coordinates yet.");
325 ebf_faces.size() < UINT_MAX,
326 "You've created a mystical element that has more faces than can be held by unsigned " 327 "int. I applaud you.");
328 const auto num_ebfs =
static_cast<unsigned int>(ebf_faces.size());
336 const unsigned int sys_dim =
337 lm_dim + num_ebfs + lm_dim *
static_cast<unsigned int>(fdf_grad_centroid_coeffs.size());
338 DenseVector<ADReal>
x(sys_dim),
b(sys_dim);
339 DenseMatrix<ADReal>
A(sys_dim, sys_dim);
342 for (
const auto lm_dim_index :
make_range(lm_dim))
345 A(lm_dim_index, lm_dim_index) = 1;
348 for (
const auto ebf_index :
make_range(num_ebfs))
349 A(lm_dim_index, lm_dim + ebf_index) = grad_ebf_coeffs[ebf_index](lm_dim_index) /
volume;
352 b(lm_dim_index) = grad_b(lm_dim_index);
355 unsigned int num_fdf_faces = 0;
358 for (
const auto ebf_index :
make_range(num_ebfs))
361 A(lm_dim + ebf_index, lm_dim + ebf_index) = 1;
363 const bool fdf_face = ebf_faces[ebf_index].second;
364 const unsigned int starting_j_index =
365 fdf_face ? lm_dim + num_ebfs + num_fdf_faces * lm_dim : 0;
367 num_fdf_faces += fdf_face;
370 for (
const auto lm_dim_index :
make_range(lm_dim))
371 A(lm_dim + ebf_index, starting_j_index + lm_dim_index) =
372 ebf_grad_coeffs[ebf_index](lm_dim_index);
375 b(lm_dim + ebf_index) = *ebf_b[ebf_index];
378 mooseAssert(num_fdf_faces == fdf_grad_centroid_coeffs.size(),
379 "Bad math in INSFVVelocityVariable::adGradlnSln(const Elem *). Please contact a " 383 for (
const auto fdf_face_index :
make_range(num_fdf_faces))
385 const auto starting_i_index = lm_dim + num_ebfs + fdf_face_index * lm_dim;
387 for (
const auto lm_dim_i_index :
make_range(lm_dim))
389 auto i_index = starting_i_index + lm_dim_i_index;
390 A(i_index, i_index) = 1;
392 for (
const auto lm_dim_j_index :
make_range(lm_dim))
394 A(i_index, lm_dim_j_index) =
395 fdf_grad_centroid_coeffs[fdf_face_index](lm_dim_i_index, lm_dim_j_index);
400 for (
const auto lm_dim_index :
make_range(lm_dim))
401 grad(lm_dim_index) =
x(lm_dim_index);
407 mooseAssert(pr.second,
"Insertion should have just happened.");
408 return pr.first->second;
417 "I believe we should only get singular systems when two-term boundary expansion is " bool _two_term_boundary_expansion
void loopOverElemFaceInfo(const Elem &elem, const MooseMesh &mesh, ActionFunctor &act, const Moose::CoordinateSystemType coord_type, const unsigned int rz_radial_coord=libMesh::invalid_uint)
unsigned int getAxisymmetricRadialCoord() const
registerMooseObject("NavierStokesApp", INSFVVelocityVariable)
const bool _cache_cell_gradients
const unsigned int invalid_uint
std::tuple< const Elem *, const Elem *, bool > determineElemOneAndTwo(const FaceInfo &fi, const FVVar &var)
const Elem & elem() const
INSFVVelocityVariable(const InputParameters ¶ms)
const Point & faceCentroid() const
void coordTransformFactor(const SubProblem &s, SubdomainID sub_id, const P &point, C &factor, SubdomainID neighbor_sub_id=libMesh::Elem::invalid_subdomain_id)
const Point & neighborCentroid() const
VectorValue< ADReal > _temp_cell_gradient
DualNumber< Real, DNDerivativeType, true > ADReal
ADReal getElemValue(const Elem *elem, const StateArg &state) const
bool isExtrapolatedBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &time) const override
virtual bool isDirichletBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
bool isInternalFace(const FaceInfo &) const
static InputParameters validParams()
void libmesh_ignore(const Args &...)
const Point & elemCentroid() const
virtual VectorValue< ADReal > uncorrectedAdGradSln(const FaceInfo &fi, const StateArg &state, const bool correct_skewness=false) const
const std::vector< double > x
VectorValue< ADReal > uncorrectedAdGradSln(const FaceInfo &fi, const StateArg &time, const bool correct_skewness=false) const override
const Point & normal() const
virtual ADReal getExtrapolatedBoundaryFaceValue(const FaceInfo &fi, bool two_term_expansion, bool correct_skewness, const Elem *elem_side_to_extrapolate_from, const StateArg &state) const
std::string grad(const std::string &var)
static InputParameters validParams()
ADReal getExtrapolatedBoundaryFaceValue(const FaceInfo &fi, bool two_term_expansion, bool correct_skewness, const Elem *elem_side_to_extrapolate_from, const StateArg &time) const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
virtual ADReal getDirichletBoundaryFaceValue(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
Moose::CoordinateSystemType getCoordSystem(SubdomainID sid) const
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
bool isFullyDevelopedFlowFace(const FaceInfo &fi) const
Returns whether the passed-in FaceInfo corresponds to a fully-developed flow face.
std::unordered_map< const Elem *, VectorValue< ADReal > > _elem_to_grad
const ADTemplateVariableGradient< Real > & adGradSln() const override