kineticGasEvaporation.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) 2017-2021 OpenCFD Ltd.
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 "kineticGasEvaporation.H"
29 #include "constants.H"
30 #include "cutCellIso.H"
31 #include "volPointInterpolation.H"
32 #include "wallPolyPatch.H"
33 #include "fvcSmooth.H"
34 
35 using namespace Foam::constant;
36 
37 
38 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
39 
40 template<class Thermo, class OtherThermo>
43 {
44  const fvMesh& mesh = this->mesh_;
45 
46  const volScalarField& alpha = this->pair().from();
47 
48  scalarField ap
49  (
51  );
52 
53  cutCellIso cutCell(mesh, ap);
54 
55  forAll(interfaceArea_, celli)
56  {
57  label status = cutCell.calcSubCell(celli, isoAlpha_);
58  interfaceArea_[celli] = 0;
59  if (status == 0) // cell is cut
60  {
61  interfaceArea_[celli] =
62  mag(cutCell.faceArea())/mesh.V()[celli];
63  }
64  }
65 
66  for (const polyPatch& pp : mesh.boundaryMesh())
67  {
68  if (isA<wallPolyPatch>(pp))
69  {
70  forAll(pp.faceCells(), faceI)
71  {
72  const label pCelli = pp.faceCells()[faceI];
73  bool interface(false);
74  if
75  (
76  sign(C_.value()) > 0
77  && (T[pCelli] - Tactivate_.value()) > 0
78  )
79  {
80  interface = true;
81  }
82 
83  if
84  (
85  sign(C_.value()) < 0
86  && (T[pCelli] - Tactivate_.value()) < 0
87  )
88  {
89  interface = true;
90  }
91 
92  if
93  (
94  interface
95  &&
96  (
97  alpha[pCelli] < 2*isoAlpha_
98  && alpha[pCelli] > 0.5*isoAlpha_
99  )
100  )
101  {
102  interfaceArea_[pCelli] =
103  mag(pp.faceAreas()[faceI])/mesh.V()[pCelli];
104  }
105  }
106  }
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
113 template<class Thermo, class OtherThermo>
116 (
117  const dictionary& dict,
118  const phasePair& pair
119 )
120 :
122  C_("C", dimless, dict),
123  Tactivate_("Tactivate", dimTemperature, dict),
124  Mv_("Mv", dimMass/dimMoles, -1, dict),
125  interfaceArea_
126  (
127  IOobject
128  (
129  "interfaceArea",
130  this->mesh_.time().timeName(),
131  this->mesh_,
134  ),
135  this->mesh_,
137  ),
138  htc_
139  (
140  IOobject
141  (
142  "htc",
143  this->mesh_.time().timeName(),
144  this->mesh_,
147  ),
148  this->mesh_,
150  ),
151  mDotc_
152  (
153  IOobject
154  (
155  "mDotc",
156  this->mesh_.time().timeName(),
157  this->mesh_,
160  ),
161  this->mesh_,
163  ),
164  isoAlpha_(dict.getOrDefault<scalar>("isoAlpha", 0.5))
165 {
166  word speciesName = IOobject::member(this->transferSpecie());
167 
168  // Get the "to" thermo
169  const typename OtherThermo::thermoType& toThermo =
170  this->getLocalThermo
171  (
172  speciesName,
173  this->toThermo_
174  );
175 
176  // Convert from g/mol to Kg/mol
177  Mv_.value() = toThermo.W()*1e-3;
178 
179  if (Mv_.value() == -1)
180  {
182  << " Please provide the molar weight (Mv) of vapour [g/mol] "
183  << abort(FatalError);
184  }
185 }
186 
187 
188 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
189 
190 template<class Thermo, class OtherThermo>
194 {
195 
196  const fvMesh& mesh = this->mesh_;
197 
198  const dimensionedScalar HerztKnudsConst
199  (
200  sqrt
201  (
202  2.0*mathematical::pi
203  * pow3(Tactivate_)
205  )
206  );
207 
208  word speciesName = IOobject::member(this->transferSpecie());
209  tmp<volScalarField> L = this->L(speciesName, T);
210 
211  updateInterface(T);
212 
213  tmp<volScalarField> tRhov
214  (
215  new volScalarField
216  (
217  IOobject
218  (
219  "tRhov",
220  mesh.time().timeName(),
221  mesh
222  ),
223  mesh,
225  )
226  );
227  volScalarField& rhov = tRhov.ref();
228 
229  tmp<volScalarField> tdeltaT
230  (
231  new volScalarField
232  (
233  IOobject
234  (
235  "tdeltaT",
236  mesh.time().timeName(),
237  mesh
238  ),
239  mesh,
241  )
242  );
243  volScalarField& deltaT = tdeltaT.ref();
244 
246 
247  if (sign(C_.value()) > 0)
248  {
249  rhov = this->pair().to().rho();
250  deltaT = max(T - Tactivate_, T0);
251  }
252  else
253  {
254  rhov = this->pair().from().rho();
255  deltaT = max(Tactivate_ - T, T0);
256  }
257 
258  htc_ = 2*mag(C_)/(2-mag(C_))*(L()*rhov/HerztKnudsConst);
259 
260  mDotc_ = htc_*deltaT*interfaceArea_;
261 
262  return tmp<volScalarField>(new volScalarField(mDotc_));
263 }
264 
265 
266 template<class Thermo, class OtherThermo>
269 (
270  label variable,
271  const volScalarField& refValue
272 )
273 {
274  if (this->modelVariable_ == variable)
275  {
276  const volScalarField coeff(htc_*interfaceArea_);
277 
278  if (sign(C_.value()) > 0)
279  {
280  return(coeff*pos(refValue - Tactivate_));
281  }
282  else
283  {
284  return(coeff*pos(Tactivate_ - refValue));
285  }
286  }
287  else
288  {
289  return tmp<volScalarField> ();
290  }
291 }
292 
293 
294 template<class Thermo, class OtherThermo>
297 (
298  label variable,
299  const volScalarField& refValue
300 )
301 {
302  if (this->modelVariable_ == variable)
303  {
304  const volScalarField coeff(htc_*interfaceArea_*Tactivate_);
305 
306  if (sign(C_.value()) > 0)
307  {
308  return(-coeff*pos(refValue - Tactivate_));
309  }
310  else
311  {
312  return(coeff*pos(Tactivate_ - refValue));
313  }
314  }
315  else
316  {
317  return tmp<volScalarField> ();
318  }
319 }
320 
321 
322 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
L
const vector L(dict.get< vector >("L"))
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:53
Foam::InterfaceCompositionModel
Base class for interface composition models, templated on the two thermodynamic models either side of...
Definition: InterfaceCompositionModel.H:58
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::meltingEvaporationModels::kineticGasEvaporation::KSu
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
Definition: kineticGasEvaporation.C:297
Foam::meltingEvaporationModels::kineticGasEvaporation
Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kin...
Definition: kineticGasEvaporation.H:205
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimDensity
const dimensionSet dimDensity
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
interface
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Foam::meltingEvaporationModels::kineticGasEvaporation::kineticGasEvaporation
kineticGasEvaporation(const dictionary &dict, const phasePair &pair)
Construct from components.
Definition: kineticGasEvaporation.C:116
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::New
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::phasePair::to
virtual const phaseModel & to() const
To phase.
Definition: phasePair.C:59
wallPolyPatch.H
Foam::dimMoles
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::constant::physicoChemical::R
const dimensionedScalar R
Universal gas constant: default SI units: [J/mol/K].
Definition: evaluateNearWall.H:6
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
cutCellIso.H
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
constants.H
Foam::meltingEvaporationModels::kineticGasEvaporation::KSp
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
Definition: kineticGasEvaporation.C:269
volPointInterpolation.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
kineticGasEvaporation.H
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::IOobject::member
word member() const
Return member (name without the extension)
Definition: IOobjectI.H:77
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::phasePair::from
virtual const phaseModel & from() const
From phase.
Definition: phasePair.C:49
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::phaseModel::rho
const dimensionedScalar & rho() const
Definition: phaseModel.H:167
Foam::meltingEvaporationModels::kineticGasEvaporation::Kexp
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
Definition: kineticGasEvaporation.C:193
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
T0
scalar T0
Definition: createFields.H:22
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177