diffusionGasEvaporation.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) 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 
29 #include "constants.H"
30 #include "cutCellIso.H"
31 #include "volPointInterpolation.H"
32 #include "fvcGrad.H"
33 
34 using namespace Foam::constant;
35 
36 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
37 
38 template<class Thermo, class OtherThermo>
41 (
42  const volScalarField& T
43 )
44 {
45  const fvMesh& mesh = this->mesh_;
46 
47  const volScalarField& alpha = this->pair().from();
48 
49  scalarField ap
50  (
52  );
53 
54  cutCellIso cutCell(mesh, ap);
55 
56  forAll(interfaceArea_, celli)
57  {
58  const label status = cutCell.calcSubCell(celli, isoAlpha_);
59  interfaceArea_[celli] = 0;
60  if (status == 0) // cell is cut
61  {
62  interfaceArea_[celli] = 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  saturationModelPtr_
80  (
82  (
83  dict.subDict("saturationPressure"),
84  this->mesh_
85  )
86  ),
87  isoAlpha_(dict.getOrDefault<scalar>("isoAlpha", 0.5)),
88  C_("C", dimless, dict),
89  Tactivate_("Tactivate", dimTemperature, 0, dict),
90  interfaceArea_
91  (
92  IOobject
93  (
94  "interfaceArea",
95  this->mesh_.time().timeName(),
96  this->mesh_,
97  IOobject::NO_READ,
98  IOobject::NO_WRITE
99  ),
100  this->mesh_,
102  ),
103  mDotc_
104  (
105  IOobject
106  (
107  "mDotc",
108  this->mesh_.time().timeName(),
109  this->mesh_,
110  IOobject::NO_READ,
111  IOobject::AUTO_WRITE
112  ),
113  this->mesh_,
115  )
116 {}
117 
119 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
120 
121 template<class Thermo, class OtherThermo>
124 ::Kexp
125 (
126  const volScalarField& T
127 )
128 {
129  const fvMesh& mesh = this->mesh_;
130 
131  const word speciesName(IOobject::member(this->transferSpecie()));
132 
133  const typename OtherThermo::thermoType& vapourThermo =
134  this->getLocalThermo
135  (
136  speciesName,
137  this->toThermo_
138  );
139 
140  const volScalarField& from = this->pair().from();
141  const volScalarField& to = this->pair().to();
142 
143  const volScalarField& Yv = this->toThermo_.composition().Y(speciesName);
144 
145  updateInterface(T);
146 
147  auto tRhog = tmp<volScalarField>::New
148  (
149  IOobject
150  (
151  "tRhog",
152  mesh.time().timeName(),
153  mesh
154  ),
155  mesh,
157  );
158  volScalarField& rhog = tRhog.ref();
159  rhog = this->pair().to().rho();
160 
161  auto tDvg = tmp<volScalarField>::New
162  (
163  IOobject
164  (
165  "tDvg",
166  mesh.time().timeName(),
167  mesh
168  ),
169  mesh,
171  );
172  volScalarField& Dvg = tDvg.ref();
173  Dvg = this->Dto(speciesName);
174 
175  tmp<volScalarField> tpSat = saturationModelPtr_->pSat(T);
176 
177  const volScalarField XvSat(tpSat()/this->toThermo_.p());
178 
179  const dimensionedScalar Wv("Wv", dimMass/dimMoles, vapourThermo.W());
180  const volScalarField YvSat
181  (
182  XvSat
183  *(
184  Wv/(XvSat*Wv + (1-XvSat)*this->toThermo_.W())
185  )
186  );
187 
188  const volScalarField Ygm(max(from*YvSat + to*Yv, Zero));
189 
190  const multiphaseInterSystem& fluid = this->fluid();
191 
192  tmp<volVectorField> tnHatInt(fluid.nVolHatfv(to, from));
193 
194  const volScalarField gradYgm(fvc::grad(Ygm) & tnHatInt());
195 
196  mDotc_ =
197  -pos(T - Tactivate_)
198  *C_*rhog*Dvg*gradYgm*interfaceArea_
199  /(1 - YvSat);
200 
201  if (debug && mesh.time().writeTime())
202  {
203  volScalarField pSat("pSat", saturationModelPtr_->pSat(T));
204  pSat.write();
205 
206  volScalarField YvSat1("YvSat", YvSat);
207  YvSat1.write();
208 
209  volScalarField YgmDebug("Ygm", Ygm);
210  YgmDebug.write();
211 
212  volScalarField gradYgmD("gradYgm", gradYgm);
213  gradYgmD.write();
214 
215  volVectorField nHatIntD("nHatInt", tnHatInt());
216  nHatIntD.write();
217  }
218 
219  return tmp<volScalarField>::New(mDotc_);
220 }
221 
222 
223 template<class Thermo, class OtherThermo>
226 KSp
227 (
228  label variable,
229  const volScalarField& refValue
230 )
231 {
232  return nullptr;
233 }
234 
235 
236 template<class Thermo, class OtherThermo>
239 KSu
240 (
241  label variable,
242  const volScalarField& refValue
243 )
244 {
245  return nullptr;
246 }
247 
248 
249 //************************************************************************ //
Different types of constants.
dimensionedScalar rhog("rhog", dimDensity, transportProperties)
twoPhaseSystem & fluid
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)
Explicit mass transfer coefficient.
static const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new 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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Base class for interface composition models, templated on the two thermodynamic models either side of...
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
word member() const
Return member (name without the extension)
Definition: IOobjectI.H:207
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)
Implicit mass transfer coefficient.
const dimensionSet dimless
Dimensionless.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual const multiphaseInter::phaseModel & to() const
To phase.
Definition: phasePair.C:52
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos(const dimensionedScalar &ds)
virtual const multiphaseInter::phaseModel & from() const
From phase.
Definition: phasePair.C:42
dynamicFvMesh & mesh
Calculate the gradient of the given field.
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
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
diffusionGasEvaporation(const dictionary &dict, const phasePair &pair)
Construct from components.
int debug
Static debugging option.
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
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
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Gas diffusion based evaporation/condensation mass transfer model.
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:172
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
virtual tmp< volScalarField > Kexp(const volScalarField &field)
Explicit total mass transfer coefficient.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127