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-2022 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 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
70 template<class Thermo, class OtherThermo>
73 (
74  const dictionary& dict,
75  const phasePair& pair
76 )
77 :
78  InterfaceCompositionModel<Thermo, OtherThermo>(dict, pair),
79  C_("C", dimless, dict),
80  Tactivate_("Tactivate", dimTemperature, dict),
81  Mv_("Mv", dimMass/dimMoles, -1, dict),
82  interfaceArea_
83  (
84  IOobject
85  (
86  "interfaceArea",
87  this->mesh_.time().timeName(),
88  this->mesh_,
89  IOobject::NO_READ,
90  IOobject::NO_WRITE
91  ),
92  this->mesh_,
94  ),
95  htc_
96  (
97  IOobject
98  (
99  "htc",
100  this->mesh_.time().timeName(),
101  this->mesh_,
102  IOobject::NO_READ,
103  IOobject::NO_WRITE
104  ),
105  this->mesh_,
107  ),
108  mDotc_
109  (
110  IOobject
111  (
112  "mDotc",
113  this->mesh_.time().timeName(),
114  this->mesh_,
115  IOobject::NO_READ,
116  IOobject::AUTO_WRITE
117  ),
118  this->mesh_,
120  ),
121  isoAlpha_(dict.getOrDefault<scalar>("isoAlpha", 0.5))
122 {
123  word speciesName = IOobject::member(this->transferSpecie());
124 
125  // Get the "to" thermo
126  const typename OtherThermo::thermoType& toThermo =
127  this->getLocalThermo
128  (
129  speciesName,
130  this->toThermo_
131  );
132 
133  // Convert from g/mol to Kg/mol
134  Mv_.value() = toThermo.W()*1e-3;
135 
136  if (Mv_.value() == -1)
137  {
139  << " Please provide the molar weight (Mv) of vapour [g/mol] "
140  << abort(FatalError);
141  }
142 }
144 
145 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
146 
147 template<class Thermo, class OtherThermo>
150 ::Kexp(const volScalarField& T)
151 {
152 
153  const fvMesh& mesh = this->mesh_;
154 
155  const dimensionedScalar HerztKnudsConst
156  (
157  sqrt
158  (
159  2.0*mathematical::pi
160  * pow3(Tactivate_)
162  )
163  );
164 
165  word speciesName = IOobject::member(this->transferSpecie());
166  tmp<volScalarField> L = mag(this->L(speciesName, T));
167 
168  updateInterface(T);
169 
170  auto tRhov = volScalarField::New
171  (
172  "tRhov",
174  mesh,
176  );
177  auto& rhov = tRhov.ref();
178 
179  auto tdeltaT = volScalarField::New
180  (
181  "tdeltaT",
183  mesh,
185  );
186  auto& deltaT = tdeltaT.ref();
187 
189 
190  if (sign(C_.value()) > 0)
191  {
192  rhov = this->pair().to().rho();
193  deltaT = max(T - Tactivate_, T0);
194  }
195  else
196  {
197  rhov = this->pair().from().rho();
198  deltaT = max(Tactivate_ - T, T0);
199  }
200 
201  htc_ = 2*mag(C_)/(2-mag(C_))*(L()*rhov/HerztKnudsConst);
202 
203  mDotc_ = htc_*deltaT*interfaceArea_;
204 
205  return tmp<volScalarField>::New(mDotc_);
206 }
207 
208 
209 template<class Thermo, class OtherThermo>
212 (
213  label variable,
214  const volScalarField& refValue
215 )
216 {
217  if (this->modelVariable_ == variable)
218  {
219  const volScalarField coeff(htc_*interfaceArea_);
220 
221  if (sign(C_.value()) > 0)
222  {
223  return(coeff*pos(refValue - Tactivate_));
224  }
225  else
226  {
227  return(coeff*pos(Tactivate_ - refValue));
228  }
229  }
230 
231  return nullptr;
232 }
233 
234 
235 template<class Thermo, class OtherThermo>
238 (
239  label variable,
240  const volScalarField& refValue
241 )
242 {
243  if (this->modelVariable_ == variable)
244  {
245  const volScalarField coeff(htc_*interfaceArea_*Tactivate_);
246 
247  if (sign(C_.value()) > 0)
248  {
249  return(-coeff*pos(refValue - Tactivate_));
250  }
251  else
252  {
253  return(coeff*pos(Tactivate_ - refValue));
254  }
255  }
256 
257  return nullptr;
258 }
259 
260 
261 // ************************************************************************* //
Different types of constants.
dimensionedScalar sign(const dimensionedScalar &ds)
const Type & value() const noexcept
Return const reference to value.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Base class for interface composition models, templated on the two thermodynamic models either side of...
dimensionedScalar sqrt(const dimensionedScalar &ds)
word member() const
Return member (name without the extension)
Definition: IOobjectI.H:207
const dimensionedScalar R
Universal gas constant: default SI units: [J/mol/K].
const dimensionSet dimless
Dimensionless.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const OtherThermo & toThermo_
Other Thermo (to)
virtual const multiphaseInter::phaseModel & to() const
To phase.
Definition: phasePair.C:52
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos(const dimensionedScalar &ds)
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
const word transferSpecie() const
Return the transferring species name.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
virtual const multiphaseInter::phaseModel & from() const
From phase.
Definition: phasePair.C:42
dynamicFvMesh & mesh
Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kin...
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
constexpr scalar pi(M_PI)
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
tmp< volScalarField > rho() const
Return the phase density.
Definition: phaseModel.C:115
const dimensionSet dimDensity
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
kineticGasEvaporation(const dictionary &dict, const phasePair &pair)
Construct from components.
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:180
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
scalar T0
Definition: createFields.H:22
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127