HeatAndMassTransferPhaseSystem.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) 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 
27 
28 #include "BlendedInterfacialModel.H"
29 #include "heatTransferModel.H"
30 #include "massTransferModel.H"
31 
32 #include "HashPtrTable.H"
33 
34 #include "fvcDiv.H"
35 #include "fvmSup.H"
36 #include "fvMatrix.H"
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasePhaseSystem>
43 (
44  const fvMesh& mesh
45 )
46 :
47  BasePhaseSystem(mesh)
48 {
50  (
51  "heatTransfer",
52  heatTransferModels_
53  );
54 
56  (
57  "massTransfer",
58  massTransferModels_
59  );
60 
62  (
64  this->phasePairs_,
65  phasePairIter
66  )
67  {
68  const phasePair& pair(phasePairIter());
69 
70  if (pair.ordered())
71  {
72  continue;
73  }
74 
75  // Initialy assume no mass transfer
76 
77  dmdt_.insert
78  (
79  pair,
80  new volScalarField
81  (
82  IOobject
83  (
84  IOobject::groupName("dmdt", pair.name()),
85  this->mesh().time().timeName(),
86  this->mesh(),
89  ),
90  this->mesh(),
92  )
93  );
94 
95  dmdtExplicit_.insert
96  (
97  pair,
98  new volScalarField
99  (
100  IOobject
101  (
102  IOobject::groupName("dmdtExplicit", pair.name()),
103  this->mesh().time().timeName(),
104  this->mesh()
105  ),
106  this->mesh(),
108  )
109  );
110 
111  volScalarField H1(heatTransferModels_[pair][pair.first()]->K());
112  volScalarField H2(heatTransferModels_[pair][pair.second()]->K());
113 
114  Tf_.insert
115  (
116  pair,
117  new volScalarField
118  (
119  IOobject
120  (
121  IOobject::groupName("Tf", pair.name()),
122  this->mesh().time().timeName(),
123  this->mesh(),
126  ),
127  (
128  H1*pair.phase1().thermo().T()
129  + H2*pair.phase2().thermo().T()
130  )
131  /max
132  (
133  H1 + H2,
135  ),
136  zeroGradientFvPatchScalarField::typeName
137  )
138  );
139  Tf_[pair]->correctBoundaryConditions();
140  }
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
146 template<class BasePhaseSystem>
149 {}
150 
151 
152 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
153 
154 template<class BasePhaseSystem>
156 (
157  const phaseModel& phase
158 ) const
159 {
160  return true;
161 }
162 
163 
164 template<class BasePhaseSystem>
167 (
168  const phasePairKey& key
169 ) const
170 {
171  const scalar dmdtSign(Pair<word>::compare(dmdt_.find(key).key(), key));
172 
173  return dmdtSign**dmdt_[key];
174 }
175 
176 
177 template<class BasePhaseSystem>
180 (
181  const Foam::phaseModel& phase
182 ) const
183 {
184  tmp<volScalarField> tdmdt
185  (
186  new volScalarField
187  (
188  IOobject
189  (
190  IOobject::groupName("dmdt", phase.name()),
191  this->mesh_.time().timeName(),
192  this->mesh_
193  ),
194  this->mesh_,
196  )
197  );
198 
200  (
202  this->phasePairs_,
203  phasePairIter
204  )
205  {
206  const phasePair& pair(phasePairIter());
207 
208  if (pair.ordered())
209  {
210  continue;
211  }
212 
213  const phaseModel* phase1 = &pair.phase1();
214  const phaseModel* phase2 = &pair.phase2();
215 
216  forAllConstIter(phasePair, pair, iter)
217  {
218  if (phase1 == &phase)
219  {
220  tdmdt() += this->dmdt(pair);
221  }
222 
223  Swap(phase1, phase2);
224  }
225  }
226 
227  return tdmdt;
228 }
229 
230 
231 template<class BasePhaseSystem>
234 {
235  autoPtr<phaseSystem::momentumTransferTable>
237 
238  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
239 
240  // Source term due to mass trasfer
242  (
244  this->phasePairs_,
245  phasePairIter
246  )
247  {
248  const phasePair& pair(phasePairIter());
249 
250  if (pair.ordered())
251  {
252  continue;
253  }
254 
255  const volVectorField& U1(pair.phase1().U());
256  const volVectorField& U2(pair.phase2().U());
257 
258  const volScalarField dmdt(this->dmdt(pair));
259  const volScalarField dmdt21(posPart(dmdt));
260  const volScalarField dmdt12(negPart(dmdt));
261 
262  *eqns[pair.phase1().name()] += dmdt21*U2 - fvm::Sp(dmdt21, U1);
263  *eqns[pair.phase2().name()] -= dmdt12*U1 - fvm::Sp(dmdt12, U2);
264  }
265 
266  return eqnsPtr;
267 }
268 
269 
270 template<class BasePhaseSystem>
273 {
274  autoPtr<phaseSystem::heatTransferTable> eqnsPtr
275  (
277  );
278 
279  phaseSystem::heatTransferTable& eqns = eqnsPtr();
280 
281  forAll(this->phaseModels_, phasei)
282  {
283  const phaseModel& phase = this->phaseModels_[phasei];
284 
285  eqns.insert
286  (
287  phase.name(),
288  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
289  );
290  }
291 
292  // Heat transfer with the interface
294  (
295  heatTransferModelTable,
296  heatTransferModels_,
297  heatTransferModelIter
298  )
299  {
300  const phasePair& pair
301  (
302  this->phasePairs_[heatTransferModelIter.key()]
303  );
304 
305  const phaseModel* phase = &pair.phase1();
306  const phaseModel* otherPhase = &pair.phase2();
307 
308  const volScalarField& Tf(*Tf_[pair]);
309 
310  const volScalarField K1
311  (
312  heatTransferModelIter()[pair.first()]->K()
313  );
314  const volScalarField K2
315  (
316  heatTransferModelIter()[pair.second()]->K()
317  );
318  const volScalarField KEff
319  (
320  K1*K2
321  /max
322  (
323  K1 + K2,
325  )
326  );
327 
328  const volScalarField* K = &K1;
329  const volScalarField* otherK = &K2;
330 
331  forAllConstIter(phasePair, pair, iter)
332  {
333  const volScalarField& he(phase->thermo().he());
334  volScalarField Cpv(phase->thermo().Cpv());
335 
336  *eqns[phase->name()] +=
337  (*K)*(Tf - phase->thermo().T())
338  + KEff/Cpv*he - fvm::Sp(KEff/Cpv, he);
339 
340  Swap(phase, otherPhase);
341  Swap(K, otherK);
342  }
343  }
344 
345  // Source term due to mass transfer
347  (
349  this->phasePairs_,
350  phasePairIter
351  )
352  {
353  const phasePair& pair(phasePairIter());
354 
355  if (pair.ordered())
356  {
357  continue;
358  }
359 
360  const phaseModel& phase1 = pair.phase1();
361  const phaseModel& phase2 = pair.phase2();
362 
363  const volScalarField& he1(phase1.thermo().he());
364  const volScalarField& he2(phase2.thermo().he());
365 
366  const volScalarField& K1(phase1.K());
367  const volScalarField& K2(phase2.K());
368 
369  const volScalarField dmdt(this->dmdt(pair));
370  const volScalarField dmdt21(posPart(dmdt));
371  const volScalarField dmdt12(negPart(dmdt));
372  const volScalarField& Tf(*Tf_[pair]);
373 
374  *eqns[phase1.name()] +=
375  dmdt21*(phase1.thermo().he(phase1.thermo().p(), Tf))
376  - fvm::Sp(dmdt21, he1)
377  + dmdt21*(K2 - K1);
378 
379  *eqns[phase2.name()] -=
380  dmdt12*(phase2.thermo().he(phase2.thermo().p(), Tf))
381  - fvm::Sp(dmdt12, he2)
382  + dmdt12*(K1 - K2);
383  }
384 
385  return eqnsPtr;
386 }
387 
388 
389 template<class BasePhaseSystem>
391 {
392  if (BasePhaseSystem::read())
393  {
394  bool readOK = true;
395 
396  // Models ...
397 
398  return readOK;
399  }
400  else
401  {
402  return false;
403  }
404 }
405 
406 
407 // ************************************************************************* //
U2
volVectorField & U2
Definition: createFields.H:23
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::phaseSystem::phasePairTable
HashTable< autoPtr< phasePair >, phasePairKey, phasePairKey::hash > phasePairTable
Definition: phaseSystem.H:133
K1
#define K1
Definition: SHA1.C:169
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::phaseSystem::heatTransferTable
HashPtrTable< fvScalarMatrix, word, string::hash > heatTransferTable
Definition: phaseSystem.H:109
Foam::dimDensity
const dimensionSet dimDensity
Foam::HeatAndMassTransferPhaseSystem::momentumTransfer
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer() const
Return the momentum transfer matrices.
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::posPart
dimensionedScalar posPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:212
Foam::phaseModel::K
virtual const volScalarField & K() const
Return the phase kinetic energy.
Foam::phaseModel::thermo
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
Foam::HeatAndMassTransferPhaseSystem::dmdt
virtual tmp< volScalarField > dmdt(const phasePairKey &key) const
Return the interfacial mass flow rate.
fvcDiv.H
Calculate the divergence of the given field.
HeatAndMassTransferPhaseSystem.H
fvMatrix.H
he2
volScalarField & he2
Definition: EEqns.H:3
Foam::heatTransferModel::dimK
static const dimensionSet dimK
Coefficient dimensions.
Definition: heatTransferModel.H:88
Foam::HeatAndMassTransferPhaseSystem::read
virtual bool read()
Read base phaseProperties dictionary.
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
massTransferModel.H
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::phaseSystem::generatePairsAndSubModels
void generatePairsAndSubModels(const word &modelName, HashTable< autoPtr< modelType >, phasePairKey, phasePairKey::hash > &models)
Generate pairs and sub-model tables.
Foam::HeatAndMassTransferPhaseSystem::heatTransfer
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:109
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
K2
#define K2
Definition: SHA1.C:170
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
U1
volVectorField & U1
Definition: createFields.H:18
Foam::HeatAndMassTransferPhaseSystem::transfersMass
virtual bool transfersMass(const phaseModel &phase) const
Return true if there is mass transfer for phase.
fvmSup.H
Calculate the matrix for implicit and explicit sources.
momentumTransfer
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:223
he
volScalarField & he
Definition: YEEqn.H:56
phasei
label phasei
Definition: pEqn.H:37
Foam::IOobject::IOobject
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:116
Foam::phaseSystem::momentumTransferTable
HashPtrTable< fvVectorMatrix, word, string::hash > momentumTransferTable
Definition: phaseSystem.H:100
HashPtrTable.H
phase2
phaseModel & phase2
Definition: createFields.H:13
Foam::basicThermo::p
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:466
Foam::phaseSystem::phasePairs_
phasePairTable phasePairs_
Phase pairs.
Definition: phaseSystem.H:167
Foam::HeatAndMassTransferPhaseSystem::~HeatAndMassTransferPhaseSystem
virtual ~HeatAndMassTransferPhaseSystem()
Destructor.
Foam::Swap
void Swap(T &a, T &b)
Definition: Swap.H:43
Foam::basicThermo::he
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
Foam::HeatAndMassTransferPhaseSystem::HeatAndMassTransferPhaseSystem
HeatAndMassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
Foam::Pair< word >::compare
static int compare(const Pair< Type > &a, const Pair< Type > &b)
Compare Pairs.
Definition: Pair.H:141
phase1
phaseModel & phase1
Definition: createFields.H:12