randomCoalescence.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 "randomCoalescence.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
35 namespace IATEsources
36 {
37  defineTypeNameAndDebug(randomCoalescence, 0);
38  addToRunTimeSelectionTable(IATEsource, randomCoalescence, dictionary);
39 }
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const IATE& iate,
50  const dictionary& dict
51 )
52 :
53  IATEsource(iate),
54  Crc_("Crc", dimless, dict),
55  C_("C", dimless, dict),
56  alphaMax_("alphaMax", dimless, dict)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
64 {
65  tmp<volScalarField> tR
66  (
67  new volScalarField
68  (
69  IOobject
70  (
71  "R",
72  iate_.phase().U().time().timeName(),
73  iate_.phase().mesh()
74  ),
75  iate_.phase().U().mesh(),
77  )
78  );
79 
80  volScalarField R = tR();
81 
82  scalar Crc = Crc_.value();
83  scalar C = C_.value();
84  scalar alphaMax = alphaMax_.value();
85  volScalarField Ut(this->Ut());
86  const volScalarField& alpha = phase();
87  const volScalarField& kappai = iate_.kappai();
88  scalar cbrtAlphaMax = cbrt(alphaMax);
89 
90  forAll(R, celli)
91  {
92  if (alpha[celli] < alphaMax - SMALL)
93  {
94  scalar cbrtAlphaMaxMAlpha = cbrtAlphaMax - cbrt(alpha[celli]);
95 
96  R[celli] =
97  (-12)*phi()*kappai[celli]*alpha[celli]
98  *Crc
99  *Ut[celli]
100  *(1 - exp(-C*cbrt(alpha[celli]*alphaMax)/cbrtAlphaMaxMAlpha))
101  /(cbrtAlphaMax*cbrtAlphaMaxMAlpha);
102  }
103  }
104 
105  return tR;
106 }
107 
108 
109 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const phaseModel & phase() const
Definition: IATEsource.H:138
virtual tmp< volScalarField > R() const
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar exp(const dimensionedScalar &ds)
const Time & time() const
Return time.
Definition: IOobject.C:227
dimensionedScalar cbrt(const dimensionedScalar &ds)
const volVectorField & U() const
Definition: phaseModel.H:170
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const IATE & iate_
Reference to the IATE this source applies to.
Definition: IATEsource.H:61
const volScalarField & kappai() const
Return the interfacial curvature.
Definition: IATE.H:120
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
randomCoalescence(const IATE &iate, const dictionary &dict)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Mesh & mesh() const
Return mesh.
const phaseModel & phase() const
Return the phase.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: PtrList.H:54
Namespace for OpenFOAM.