LegendreMagnaudet.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) 2014-2023 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 "LegendreMagnaudet.H"
27 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace liftModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phaseInterface& interface
48 )
49 :
50  dispersedLiftModel(dict, interface),
51  residualRe_("residualRe", dimless, dict)
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
64 {
65  const volScalarField Re(max(interface_.Re(), residualRe_));
66 
67  const volScalarField Sr
68  (
69  sqr(interface_.dispersed().d())
70  /(Re*interface_.continuous().thermo().nu())
71  *mag(fvc::grad(interface_.continuous().U()))
72  );
73 
74  const volScalarField ClLowSqr
75  (
76  sqr(6*2.255)*sqr(Sr)
78  );
79 
80  const volScalarField ClHighSqr(sqr(0.5*(Re + 16)/(Re + 29)));
81 
82  return sqrt(ClLowSqr + ClHighSqr);
83 }
84 
85 
86 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Model for the lift force between two phases.
Definition: liftModel.H:53
Lift model of Legendre and Magnaudet.
LegendreMagnaudet(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
virtual tmp< volScalarField > Cl() const
Lift coefficient.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the gradient of the given field.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
addToRunTimeSelectionTable(liftModel, constantLiftCoefficient, dictionary)
defineTypeNameAndDebug(constantLiftCoefficient, 0)
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict