SchillerNaumannDrag.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "SchillerNaumannDrag.H"
29 #include "coupledToFluid.H"
30 #include "sphericalCoupled.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace Lagrangian
37 {
40  (
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 Foam::Lagrangian::SchillerNaumannDrag::calcD
53 (
54  const LagrangianModelRef& model,
55  const LagrangianSubMesh& subMesh
56 ) const
57 {
58  const clouds::spherical& sCloud = cloud<clouds::spherical>();
59  const clouds::sphericalCoupled& scCloud = cloud<clouds::sphericalCoupled>();
60 
61  const LagrangianSubScalarField& Re = scCloud.Re(model, subMesh);
62  const LagrangianSubScalarSubField d(sCloud.d(model, subMesh));
63 
65  <
66  clouds::coupledToIncompressibleFluid,
67  clouds::coupledToFluid
68  >();
69 
70  tmp<LagrangianSubScalarField> tmucByRhoOrMuc =
71  isCloud<clouds::coupledToIncompressibleFluid>()
72  ? (
73  cloud<clouds::coupledToIncompressibleFluid>().nuc(model, subMesh)
74  /cloud<clouds::coupledToIncompressibleFluid>().rhoByRhoc
75  )
76  : tmp<LagrangianSubScalarField>
77  (
78  cloud<clouds::coupledToFluid>().muc(model, subMesh)
79  );
80 
81  return
83  (
84  "D:" + Foam::name(subMesh.group()),
85  CdRe(Re)*(constant::mathematical::pi/8)*d*tmucByRhoOrMuc
86  );
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
93 (
94  const word& name,
95  const LagrangianMesh& mesh,
96  const dictionary& modelDict,
97  const dictionary& stateDict
98 )
99 :
100  drag(name, mesh, modelDict, stateDict)
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
108 (
110 )
111 {
112  return
113  neg(Re - 1000)*24*(1 + 0.15*pow(Re, 0.687))
114  + pos0(Re - 1000)*0.44*Re;
115 }
116 
117 
118 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Class containing Lagrangian geometry and topology.
Base class for Lagrangian models.
static tmp< LagrangianSubScalarField > CdRe(const LagrangianSubScalarField &Re)
Return the drag coefficient times Reynold's number, as a function.
SchillerNaumannDrag(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
void assertCloud() const
Generate an error if the cloud is not one of the given types.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
dimensionedScalar neg(const dimensionedScalar &ds)
LagrangianSubSubField< scalar > LagrangianSubScalarSubField
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
LagrangianSubField< scalar > LagrangianSubScalarField
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97