VakhrushevEfremov.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
28#include "VakhrushevEfremov.H"
29#include "phasePair.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace aspectRatioModels
37{
40 (
41 aspectRatioModel,
43 dictionary
44 );
45}
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const dictionary& dict,
54 const phasePair& pair
55)
56:
57 aspectRatioModel(dict, pair)
58{}
59
60
61// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62
64{}
65
66
67// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68
71{
72 volScalarField Ta(pair_.Ta());
73
74 return
75 neg(Ta - scalar(1))*scalar(1)
76 + pos0(Ta - scalar(1))*neg(Ta - scalar(39.8))
77 *pow3(0.81 + 0.206*tanh(1.6 - 2*log10(max(Ta, scalar(1)))))
78 + pos0(Ta - scalar(39.8))*0.24;
79}
80
81
82// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Aspect ratio model of Vakhrushev and Efremov.
Aspect ratio model of Vakhrushev and Efremov.
virtual tmp< volScalarField > E() const
Aspect ratio.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
dictionary dict