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