Frank.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 "Frank.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace wallLubricationModels
35 {
36  defineTypeNameAndDebug(Frank, 0);
38  (
39  wallLubricationModel,
40  Frank,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const phasePair& pair
53 )
54 :
55  wallLubricationModel(dict, pair),
56  Cwd_("Cwd", dimless, dict),
57  Cwc_("Cwc", dimless, dict),
58  p_(readScalar(dict.lookup("p")))
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 {
72  volVectorField Ur(pair_.Ur());
73 
74  const volVectorField& n(nWall());
75  const volScalarField& y(yWall());
76 
77  volScalarField Eo(pair_.Eo());
78  volScalarField yTilde(y/(Cwc_*pair_.dispersed().d()));
79 
80  return
81  (
82  pos(Eo - 1.0)*neg(Eo - 5.0)*exp(-0.933*Eo + 0.179)
83  + pos(Eo - 5.0)*neg(Eo - 33.0)*(0.00599*Eo - 0.0187)
84  + pos(Eo - 33.0)*0.179
85  )
86  *max
87  (
88  dimensionedScalar("zero", dimless/dimLength, 0.0),
89  (1.0 - yTilde)/(Cwd_*y*pow(yTilde, p_ - 1.0))
90  )
91  *pair_.continuous().rho()
92  *magSqr(Ur - (Ur & n)*n)
93  *n;
94 }
95 
96 
97 // ************************************************************************* //
Frank(const dictionary &dict, const phasePair &pair)
Construct from components.
dictionary dict
const volScalarField & yWall() const
tmp< volVectorField > Ur() const
Relative velocity.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const phasePair & pair_
Phase pair.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual const phaseModel & continuous() const
Continuous phase.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar pos(const dimensionedScalar &ds)
scalar y
stressControl lookup("compactNormalStress") >> compactNormalStress
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
const volVectorField & nWall() const
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > d() const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< volVectorField > Fi() const
Return phase-intensive wall lubrication force.
label n
tmp< volScalarField > Eo() const
Eotvos number.
A class for managing temporary objects.
Definition: PtrList.H:54
virtual ~Frank()
Destructor.
Namespace for OpenFOAM.