LopezDeBertodano.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-2018 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 "LopezDeBertodano.H"
27 #include "phasePair.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace turbulentDispersionModels
36 {
37  defineTypeNameAndDebug(LopezDeBertodano, 0);
39  (
40  turbulentDispersionModel,
41  LopezDeBertodano,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const dictionary& dict,
53  const phasePair& pair
54 )
55 :
56  turbulentDispersionModel(dict, pair),
57  Ctd_("Ctd", dimless, dict)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 {
72  return
73  Ctd_
74  *pair_.continuous().rho()
75  *pair_.continuous().turbulence().k();
76 }
77 
78 
79 // ************************************************************************* //
dictionary dict
const PhaseCompressibleTurbulenceModel< phaseModel > & turbulence() const
Return the turbulence model.
LopezDeBertodano(const dictionary &dict, const phasePair &pair)
Construct from a dictionary and a phase pair.
virtual tmp< volScalarField > D() const
Turbulent diffusivity.
Macros for easy insertion into run-time selection tables.
virtual const phaseModel & continuous() const
Continuous phase.
const phasePair & pair_
Phase pair.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.