https://mooseframework.inl.gov
IntegralDirectedSurfaceForce.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 
11 #include "MathFVUtils.h"
12 #include "NSFVUtils.h"
13 #include "NS.h"
14 #include "MooseFunctorArguments.h"
15 
16 #include <cmath>
17 
19 
22 {
24  params.addClassDescription(
25  "Computes the directed force coming from friction and pressure differences on a surface. One "
26  "can use this object for the computation of the drag and lift coefficient as well.");
27  params.addRequiredParam<MooseFunctorName>("vel_x", "The velocity in direction x.");
28  params.addParam<MooseFunctorName>("vel_y", "The velocity in direction y.");
29  params.addParam<MooseFunctorName>("vel_z", "The velocity in direction z.");
30  params.addRequiredParam<MooseFunctorName>(NS::mu, "The dynamic viscosity.");
31  params.addRequiredParam<MooseFunctorName>(NS::pressure, "The pressure functor.");
32  params.addRequiredParam<RealVectorValue>("principal_direction",
33  "The direction in which the force is computed.");
34  return params;
35 }
36 
38  : SideIntegralPostprocessor(parameters),
39  _mu(getFunctor<Real>(NS::mu)),
40  _pressure(getFunctor<Real>(NS::pressure)),
41  _direction(getParam<RealVectorValue>("principal_direction"))
42 {
43  _vel_components.push_back(&getFunctor<Real>("vel_x"));
44 
45  if (_mesh.dimension() >= 2)
46  {
47  if (!isParamValid("vel_y"))
48  paramError("vel_y",
49  "For 2D meshes the second velocity component should be provided as well!");
50  _vel_components.push_back(&getFunctor<Real>("vel_y"));
51 
52  if (_mesh.dimension() == 3)
53  {
54  if (!isParamValid("vel_z"))
55  paramError("vel_z",
56  "For 3D meshes the third velocity component should be provided as well!");
57  _vel_components.push_back(&getFunctor<Real>("vel_z"));
58  }
59  }
60 
61  // This object will only work with finite volume variables
62  _qp_integration = false;
63 
64  checkFunctorSupportsSideIntegration<Real>("vel_x", _qp_integration);
65  if (_vel_components.size() >= 2)
66  checkFunctorSupportsSideIntegration<Real>("vel_y", _qp_integration);
67  if (_vel_components.size() == 3)
68  checkFunctorSupportsSideIntegration<Real>("vel_z", _qp_integration);
69 }
70 
71 Real
73 {
74  mooseAssert(fi, "We should have a face info in " + name());
75 
76  const auto state = determineState();
77  const auto face_arg = Moose::FaceArg(
78  {fi, Moose::FV::LimiterType::CentralDifference, true, false, nullptr, nullptr});
79  const auto elem_arg = Moose::ElemArg({fi->elemPtr(), false});
80 
81  RealTensorValue pressure_term;
82  RealVectorValue cell_velocity = 0;
83  RealVectorValue face_velocity = 0;
84  const Real pressure = _pressure(face_arg, state);
85  const Real mu = _mu(face_arg, state);
86  for (const auto i : make_range(_mesh.dimension()))
87  {
88  cell_velocity(i) = (*_vel_components[i])(elem_arg, state);
89  face_velocity(i) = (*_vel_components[i])(face_arg, state);
90  pressure_term(i, i) = -pressure;
91  }
92 
93  const auto shear_force = mu *
94  (cell_velocity - face_velocity -
95  (cell_velocity - face_velocity) * fi->normal() * fi->normal()) /
96  std::abs(fi->dCN() * fi->normal());
97 
98  return (shear_force * _direction - pressure_term * fi->normal() * _direction);
99 }
100 
101 Real
103 {
104  mooseError(this->type() + " does not have an implementation for quadrature-based evaluation!");
105  return 0.0;
106 }
const Moose::Functor< Real > & _pressure
Pressure field.
const Moose::Functor< Real > & _mu
The dynamic viscosity.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
Moose::StateArg determineState() const
Real computeFaceInfoIntegral(const FaceInfo *fi) override
virtual const std::string & name() const
registerMooseObject("NavierStokesApp", IntegralDirectedSurfaceForce)
void addRequiredParam(const std::string &name, const std::string &doc_string)
bool isParamValid(const std::string &name) const
TensorValue< Real > RealTensorValue
virtual unsigned int dimension() const
static const std::string mu
Definition: NS.h:123
const std::string & type() const
const Point & normal() const
void paramError(const std::string &param, Args... args) const
IntegralDirectedSurfaceForce(const InputParameters &parameters)
const RealVectorValue _direction
The direction in which the force is measured.
const Elem * elemPtr() const
std::vector< const Moose::Functor< Real > * > _vel_components
Velocity components.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Postprocessor which computes the directed force coming from friction and pressure differences on a su...
static const std::string pressure
Definition: NS.h:56
IntRange< T > make_range(T beg, T end)
void mooseError(Args &&... args) const
void addClassDescription(const std::string &doc_string)
const Point & dCN() const