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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2017 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 "Frank.H"
29 #include "phasePair.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace wallLubricationModels
37 {
38  defineTypeNameAndDebug(Frank, 0);
40  (
41  wallLubricationModel,
42  Frank,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const dictionary& dict,
54  const phasePair& pair
55 )
56 :
57  wallLubricationModel(dict, pair),
58  Cwd_("Cwd", dimless, dict),
59  Cwc_("Cwc", dimless, dict),
60  p_(dict.get<scalar>("p"))
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71 
73 {
74  volVectorField Ur(pair_.Ur());
75 
76  const volVectorField& n(nWall());
77  const volScalarField& y(yWall());
78 
79  volScalarField Eo(pair_.Eo());
80  volScalarField yTilde(y/(Cwc_*pair_.dispersed().d()));
81 
82  return
83  (
84  pos0(Eo - 1.0)*neg(Eo - 5.0)*exp(-0.933*Eo + 0.179)
85  + pos0(Eo - 5.0)*neg(Eo - 33.0)*(0.00599*Eo - 0.0187)
86  + pos0(Eo - 33.0)*0.179
87  )
88  *max
89  (
91  (1.0 - yTilde)/(Cwd_*y*pow(yTilde, p_ - 1.0))
92  )
93  *pair_.continuous().rho()
94  *magSqr(Ur - (Ur & n)*n)
95  *n;
96 }
97 
98 
99 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Frank(const dictionary &dict, const phasePair &pair)
Construct from components.
Definition: Frank.C:46
dictionary dict
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
addToRunTimeSelectionTable(wallLubricationModel, Antal, dictionary)
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
scalar y
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< volVectorField > Fi() const
Return phase-intensive wall lubrication force.
Definition: Frank.C:66
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label n
virtual ~Frank()
Destructor.
Definition: Frank.C:60
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133