LegendreMagnaudet.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "LegendreMagnaudet.H"
29 #include "phasePair.H"
30 #include "fvcGrad.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace liftModels
38 {
39  defineTypeNameAndDebug(LegendreMagnaudet, 0);
40  addToRunTimeSelectionTable(liftModel, LegendreMagnaudet, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const dictionary& dict,
50  const phasePair& pair
51 )
52 :
53  liftModel(dict, pair),
54  residualRe_("residualRe", dimless, dict)
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 {
68  volScalarField Re(max(pair_.Re(), residualRe_));
69 
71  (
72  sqr(pair_.dispersed().d())
73  /(
74  Re
75  *pair_.continuous().nu()
76  )
77  *mag(fvc::grad(pair_.continuous().U()))
78  );
79 
80  volScalarField ClLowSqr
81  (
82  sqr(6*2.255)
83  *sqr(Sr)
84  /(
86  *Re
87  *pow3(Sr + 0.2*Re)
88  )
89  );
90 
91  volScalarField ClHighSqr
92  (
93  sqr(0.5*(Re + 16)/(Re + 29))
94  );
95 
96  return sqrt(ClLowSqr + ClHighSqr);
97 }
98 
99 
100 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
addToRunTimeSelectionTable(liftModel, constantLiftCoefficient, dictionary)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
defineTypeNameAndDebug(constantLiftCoefficient, 0)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
virtual ~LegendreMagnaudet()
Destructor.
virtual tmp< volScalarField > Cl() const
Lift coefficient.
Calculate the gradient of the given field.
constexpr scalar pi(M_PI)
scalarField Re(const UList< complex > &cmplx)
Extract real component.
Definition: complexField.C:207
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
Namespace for OpenFOAM.
LegendreMagnaudet(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.