turbulentBreakUp.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) 2013-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 "turbulentBreakUp.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
35 namespace IATEsources
36 {
37  defineTypeNameAndDebug(turbulentBreakUp, 0);
38  addToRunTimeSelectionTable(IATEsource, turbulentBreakUp, dictionary);
39 }
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const IATE& iate,
50  const dictionary& dict
51 )
52 :
53  IATEsource(iate),
54  Cti_("Cti", dimless, dict),
55  WeCr_("WeCr", dimless, dict)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
60 
63 {
64  tmp<volScalarField> tR
65  (
66  new volScalarField
67  (
68  IOobject
69  (
70  "R",
71  iate_.phase().time().timeName(),
72  iate_.phase().mesh()
73  ),
74  iate_.phase().mesh(),
76  )
77  );
78 
79  volScalarField R = tR();
80 
81  scalar Cti = Cti_.value();
82  scalar WeCr = WeCr_.value();
83  volScalarField Ut(this->Ut());
84  volScalarField We(this->We());
85  const volScalarField& d(iate_.d()());
86 
87  forAll(R, celli)
88  {
89  if (We[celli] > WeCr)
90  {
91  R[celli] =
92  (1.0/3.0)
93  *Cti/d[celli]
94  *Ut[celli]
95  *sqrt(1 - WeCr/We[celli])
96  *exp(-WeCr/We[celli]);
97  }
98  }
99 
100  return tR;
101 }
102 
103 
104 // ************************************************************************* //
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
R
#define R(A, B, C, D, E, F, K, M)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::diameterModels::IATEsources::turbulentBreakUp::R
virtual tmp< volScalarField > R() const
dict
dictionary dict
Definition: searchingEngine.H:14
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::diameterModels::IATEsources::turbulentBreakUp::turbulentBreakUp
turbulentBreakUp(const IATE &iate, const dictionary &dict)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)