IATE.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 "IATE.H"
27 #include "IATEsource.H"
28 #include "fvmDdt.H"
29 #include "fvmDiv.H"
30 #include "fvmSup.H"
31 #include "fvcDdt.H"
32 #include "fvcDiv.H"
33 #include "fvcAverage.H"
34 #include "fvOptions.H"
35 #include "mathematicalConstants.H"
36 #include "fundamentalConstants.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace diameterModels
44 {
45  defineTypeNameAndDebug(IATE, 0);
46 
48  (
49  diameterModel,
50  IATE,
51  dictionary
52  );
53 }
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const dictionary& diameterProperties,
62  const phaseModel& phase
63 )
64 :
65  diameterModel(diameterProperties, phase),
66  kappai_
67  (
68  IOobject
69  (
70  IOobject::groupName("kappai", phase.name()),
71  phase_.time().timeName(),
72  phase_.mesh(),
73  IOobject::MUST_READ,
74  IOobject::AUTO_WRITE
75  ),
76  phase_.mesh()
77  ),
78  dMax_("dMax", dimLength, diameterProperties_),
79  dMin_("dMin", dimLength, diameterProperties_),
80  residualAlpha_
81  (
82  "residualAlpha",
83  dimless,
84  diameterProperties_
85  ),
86  d_
87  (
88  IOobject
89  (
90  IOobject::groupName("d", phase.name()),
91  phase_.time().timeName(),
92  phase_.mesh(),
93  IOobject::NO_READ,
94  IOobject::AUTO_WRITE
95  ),
96  dsm()
97  ),
98  sources_
99  (
100  diameterProperties_.lookup("sources"),
101  IATEsource::iNew(*this)
102  )
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
107 
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::dsm() const
115 {
116  return max(6/max(kappai_, 6/dMax_), dMin_);
117 }
118 
119 // Placeholder for the nucleation/condensation model
120 // Foam::tmp<Foam::volScalarField> Foam::diameterModels::IATE::Rph() const
121 // {
122 // const volScalarField& T = phase_thermo().T();
123 // const volScalarField& p = phase_.p();
124 //
125 // scalar A, B, C, sigma, vm, Rph;
126 //
127 // volScalarField ps(1e5*pow(10, A - B/(T + C)));
128 // volScalarField Dbc
129 // (
130 // 4*sigma*vm/(constant::physicoChemical::k*T*log(p/ps))
131 // );
132 //
133 // return constant::mathematical::pi*sqr(Dbc)*Rph;
134 // }
135 
137 {
138  // Initialise the accumulated source term to the dilatation effect
140  (
141  (
142  (1.0/3.0)
143  /max
144  (
145  fvc::average(phase_ + phase_.oldTime()),
146  residualAlpha_
147  )
148  )
149  *(fvc::ddt(phase_) + fvc::div(phase_.alphaPhi()))
150  );
151 
152  // Accumulate the run-time selectable sources
153  forAll(sources_, j)
154  {
155  R -= sources_[j].R();
156  }
157 
158  fv::options& fvOptions(fv::options::New(phase_.mesh()));
159 
160  // Construct the interfacial curvature equation
161  fvScalarMatrix kappaiEqn
162  (
163  fvm::ddt(kappai_) + fvm::div(phase_.phi(), kappai_)
164  - fvm::Sp(fvc::div(phase_.phi()), kappai_)
165  ==
166  - fvm::SuSp(R, kappai_)
167  //+ Rph() // Omit the nucleation/condensation term
168  + fvOptions(kappai_)
169  );
170 
171  kappaiEqn.relax();
172 
173  fvOptions.constrain(kappaiEqn);
174 
175  kappaiEqn.solve();
176 
177  // Update the Sauter-mean diameter
178  d_ = dsm();
179 }
180 
181 
182 bool Foam::diameterModels::IATE::read(const dictionary& phaseProperties)
183 {
184  diameterModel::read(phaseProperties);
185 
186  diameterProperties_.lookup("dMax") >> dMax_;
187  diameterProperties_.lookup("dMin") >> dMin_;
188 
189  // Re-create all the sources updating number, type and coefficients
190  PtrList<IATEsource>
191  (
192  diameterProperties_.lookup("sources"),
193  IATEsource::iNew(*this)
194  ).transfer(sources_);
195 
196  return true;
197 }
198 
199 
200 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual bool read(const dictionary &phaseProperties)=0
Read phaseProperties dictionary.
void transfer(dictionary &)
Transfer the contents of the argument and annul the argument.
Definition: dictionary.C:1067
Fundamental dimensioned constants.
Macros for easy insertion into run-time selection tables.
fv::options & fvOptions
Area-weighted average a surfaceField creating a volField.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dictionary()
Construct top-level dictionary null.
Definition: dictionary.C:113
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculate the first temporal derivative.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
stressControl lookup("compactNormalStress") >> compactNormalStress
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
dynamicFvMesh & mesh
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual ~IATE()
Destructor.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Calulate the matrix for the first temporal derivative.
IATE(const dictionary &diameterProperties, const phaseModel &phase)
Construct from components.
word timeName
Definition: getTimeIndex.H:3
Calculate the divergence of the given field.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Calculate the matrix for the divergence of the given field and flux.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
const Mesh & mesh() const
Return mesh.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
const surfaceScalarField & alphaPhi() const
Definition: phaseModel.H:200
virtual void correct()
Correct the diameter field.
A class for managing temporary objects.
Definition: PtrList.H:54
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
const surfaceScalarField & phi() const
Definition: phaseModel.H:190
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.