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-2022 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 {
36  defineTypeNameAndDebug(LegendreMagnaudet, 0);
37  addToRunTimeSelectionTable(liftModel, LegendreMagnaudet, dictionary);
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  (
70  /(
71  Re
73  )
75  );
76 
77  const volScalarField ClLowSqr
78  (
79  sqr(6*2.255)
80  *sqr(Sr)
81  /(
83  *Re
84  *pow3(Sr + 0.2*Re)
85  )
86  );
87 
88  const volScalarField ClHighSqr
89  (
90  sqr(0.5*(Re + 16)/(Re + 29))
91  );
92 
93  return sqrt(ClLowSqr + ClHighSqr);
94 }
95 
96 
97 // ************************************************************************* //
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
LegendreMagnaudet(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
virtual tmp< volVectorField > U() const =0
Return the velocity.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dispersedPhaseInterface interface_
Interface.
const dimensionSet dimless
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > Cl() const
Lift coefficient.
Calculate the gradient of the given field.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
virtual ~LegendreMagnaudet()
Destructor.
const phaseModel & continuous() const
Continuous phase.
const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
tmp< volScalarField > Re() const
Reynolds number.
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:73
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.