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 | Copyright (C) 2014-2015 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 "phasePair.H"
28 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace liftModels
36 {
37  defineTypeNameAndDebug(LegendreMagnaudet, 0);
38  addToRunTimeSelectionTable(liftModel, LegendreMagnaudet, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
48  const phasePair& pair
49 )
50 :
51  liftModel(dict, pair),
52  residualRe_("residualRe", dimless, dict)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
65 {
66  volScalarField Re(max(pair_.Re(), residualRe_));
67 
69  (
70  sqr(pair_.dispersed().d())
71  /(
72  Re
73  *pair_.continuous().nu()
74  )
76  );
77 
78  volScalarField ClLowSqr
79  (
80  sqr(6.0*2.255)
81  *sqr(Sr)
82  /(
84  *Re
85  *pow3(Sr + 0.2*Re)
86  )
87  );
88 
89  volScalarField ClHighSqr
90  (
91  sqr(0.5*(Re + 16.0)/(Re + 29.0))
92  );
93 
94  return sqrt(ClLowSqr + ClHighSqr);
95 }
96 
97 
98 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const phasePair & pair_
Phase pair.
Definition: liftModel.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual const phaseModel & continuous() const
Continuous phase.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the gradient of the given field.
virtual ~LegendreMagnaudet()
Destructor.
virtual tmp< volScalarField > Cl() const
Lift coefficient.
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
tmp< volScalarField > Re() const
Reynolds number.
const volVectorField & U() const
Definition: phaseModel.H:170
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > d() const
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:54
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.
LegendreMagnaudet(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.