randomCoalescence.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) 2013-2019 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"
27 #include "fvmSup.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace IATEsources
37 {
38  defineTypeNameAndDebug(randomCoalescence, 0);
39  addToRunTimeSelectionTable(IATEsource, randomCoalescence, dictionary);
40 }
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
49 (
50  const IATE& iate,
51  const dictionary& dict
52 )
53 :
54  IATEsource(iate),
55  Crc_("Crc", dimless, dict),
56  C_("C", dimless, dict),
57  alphaMax_("alphaMax", dimless, dict)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
65 (
66  const volScalarField& alphai,
67  volScalarField& kappai
68 ) const
69 {
71  (
72  IOobject
73  (
74  "randomCoalescence:R",
75  iate_.phase().time().timeName(),
76  iate_.phase().mesh()
77  ),
78  iate_.phase().mesh(),
80  );
81 
82  const scalar Crc = Crc_.value();
83  const scalar C = C_.value();
84  const scalar alphaMax = alphaMax_.value();
85  const volScalarField Ut(this->Ut());
86  const volScalarField& alpha = phase();
87  const scalar cbrtAlphaMax = cbrt(alphaMax);
88 
89  forAll(R, celli)
90  {
91  if (kappai[celli] > 0 && alpha[celli] < alphaMax - small)
92  {
93  const scalar cbrtAlphaMaxMAlpha = cbrtAlphaMax - cbrt(alpha[celli]);
94 
95  R[celli] =
96  12*phi()*kappai[celli]*alpha[celli]
97  *Crc
98  *Ut[celli]
99  *(1 - exp(-C*cbrt(alpha[celli]*alphaMax)/cbrtAlphaMaxMAlpha))
100  /(cbrtAlphaMax*cbrtAlphaMaxMAlpha);
101  }
102  }
103 
104  return -fvm::Sp(R, kappai);
105 }
106 
107 
108 // ************************************************************************* //
dictionary dict
tmp< volScalarField > Ut() const
Return the bubble turbulent velocity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volScalarField > R() const
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const phaseModel & phase() const
Return the phase.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
const phaseModel & phase() const
Definition: IATEsource.H:138
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
const IATE & iate_
Reference to the IATE this source applies to.
Definition: IATEsource.H:61
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 Time & time() const
Return time.
Definition: IOobject.C:360
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
virtual tmp< volScalarField > R() const =0
Namespace for OpenFOAM.