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 | 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 "VakhrushevEfremov.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace aspectRatioModels
35 {
36  defineTypeNameAndDebug(VakhrushevEfremov, 0);
38  (
39  aspectRatioModel,
40  VakhrushevEfremov,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const phasePair& pair
53 )
54 :
55  aspectRatioModel(dict, pair)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
69 {
70  volScalarField Ta(pair_.Ta());
71 
72  return
73  neg(Ta - scalar(1))*scalar(1)
74  + pos(Ta - scalar(1))*neg(Ta - scalar(39.8))
75  *pow3(0.81 + 0.206*tanh(1.6 - 2*log10(max(Ta, scalar(1)))))
76  + pos(Ta - scalar(39.8))*0.24;
77 }
78 
79 
80 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
dictionary dict
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
const phasePair & pair_
Phase pair.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar pos(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > Ta() const
Takahashi Number.
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > E() const
Aspect ratio.
A class for managing temporary objects.
Definition: PtrList.H:54
VakhrushevEfremov(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and an ordered phase pair.
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.