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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 "randomCoalescence.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace IATEsources
38 {
39  defineTypeNameAndDebug(randomCoalescence, 0);
40  addToRunTimeSelectionTable(IATEsource, randomCoalescence, dictionary);
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
50 (
51  const IATE& iate,
52  const dictionary& dict
53 )
54 :
55  IATEsource(iate),
56  Crc_("Crc", dimless, dict),
57  C_("C", dimless, dict),
58  alphaMax_("alphaMax", dimless, dict)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
66 {
68  (
69  new volScalarField
70  (
71  IOobject
72  (
73  "R",
74  iate_.phase().U().time().timeName(),
75  iate_.phase().mesh()
76  ),
77  iate_.phase().U().mesh(),
79  )
80  );
81 
82  volScalarField R = tR();
83 
84  scalar Crc = Crc_.value();
85  scalar C = C_.value();
86  scalar alphaMax = alphaMax_.value();
87  volScalarField Ut(this->Ut());
88  const volScalarField& alpha = phase();
89  const volScalarField& kappai = iate_.kappai();
90  scalar cbrtAlphaMax = cbrt(alphaMax);
91 
92  forAll(R, celli)
93  {
94  if (alpha[celli] < alphaMax - SMALL)
95  {
96  scalar cbrtAlphaMaxMAlpha = cbrtAlphaMax - cbrt(alpha[celli]);
97 
98  R[celli] =
99  (-12)*phi()*kappai[celli]*alpha[celli]
100  *Crc
101  *Ut[celli]
102  *(1 - exp(-C*cbrt(alpha[celli]*alphaMax)/cbrtAlphaMaxMAlpha))
103  /(cbrtAlphaMax*cbrtAlphaMaxMAlpha);
104  }
105  }
106 
107  return tR;
108 }
109 
110 
111 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
Graphite solid properties.
Definition: C.H:46
const volScalarField & kappai() const
Return the interfacial curvature.
Definition: IATE.H:139
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const volVectorField & U() const
Definition: phaseModel.H:198
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
const phaseModel & phase() const
Return the phase.
addToRunTimeSelectionTable(IATEsource, dummy, word)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
const phaseModel & phase() const
Definition: IATEsource.H:145
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
Definition: IATEsource.C:84
const Mesh & mesh() const noexcept
Return mesh.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const IATE & iate_
Reference to the IATE this source applies to.
Definition: IATEsource.H:62
randomCoalescence(const IATE &iate, const dictionary &dict)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127