Frank.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-2023 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace wallLubricationModels
34 {
37  (
39  Frank,
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const phaseInterface& interface
52 )
53 :
55  Cwd_("Cwd", dimless, dict),
56  Cwc_("Cwc", dimless, dict),
57  p_(dict.lookup<scalar>("p"))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
70 {
71  const volVectorField Ur(interface_.Ur());
72 
73  const volVectorField& n(nWall());
74  const volScalarField& y(yWall());
75 
76  const volScalarField Eo(interface_.Eo());
77  const volScalarField yTilde(y/(Cwc_*interface_.dispersed().d()));
78 
79  return zeroGradWalls
80  (
81  (
82  pos0(Eo - 1)*neg(Eo - 5)*exp(-0.933*Eo + 0.179)
83  + pos0(Eo - 5)*neg(Eo - 33)*(0.00599*Eo - 0.0187)
84  + pos0(Eo - 33)*0.179
85  )
86  *max
87  (
89  (1 - yTilde)/(Cwd_*y*pow(yTilde, p_ - 1))
90  )
91  *interface_.continuous().rho()
92  *magSqr(Ur - (Ur & n)*n)
93  *n
94  );
95 }
96 
97 
98 // ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
Model for the wall lubrication force between two phases.
Wall lubrication model of Frank.
Definition: Frank.H:69
Frank(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
Definition: Frank.C:49
virtual ~Frank()
Destructor.
Definition: Frank.C:63
tmp< volVectorField > Fi() const
Return phase-intensive wall lubrication force.
Definition: Frank.C:69
addToRunTimeSelectionTable(wallLubricationModel, Antal, dictionary)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict