TomiyamaWallLubrication.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-------------------------------------------------------------------------------
10License
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
29#include "phasePair.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace wallLubricationModels
37{
38 defineTypeNameAndDebug(TomiyamaWallLubrication, 0);
40 (
41 wallLubricationModel,
42 TomiyamaWallLubrication,
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 D_("Cwd", dimLength, dict)
59{}
60
61
62// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63
65{}
66
67
68// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69
72{
73 volVectorField Ur(pair_.Ur());
74
75 const volVectorField& n(nWall());
76 const volScalarField& y(yWall());
77
78 volScalarField Eo(pair_.Eo());
79
80 return
81 (
82 pos0(Eo - 1.0)*neg(Eo - 5.0)*exp(-0.933*Eo + 0.179)
83 + pos0(Eo - 5.0)*neg(Eo - 33.0)*(0.00599*Eo - 0.0187)
84 + pos0(Eo - 33.0)*0.179
85 )
86 *0.5
87 *pair_.dispersed().d()
88 *(
89 1/sqr(y)
90 - 1/sqr(D_ - y)
91 )
92 *pair_.continuous().rho()
93 *magSqr(Ur - (Ur & n)*n)
94 *n;
95}
96
97
98// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A class for managing temporary objects.
Definition: tmp.H:65
tmp< volVectorField > Fi() const
Return phase-intensive wall lubrication force.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedScalar neg(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict