Newer
Older
Henry Weller
committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "HeatAndMassTransferPhaseSystem.H"
#include "BlendedInterfacialModel.H"
#include "heatTransferModel.H"
#include "massTransferModel.H"
#include "interfaceCompositionModel.H"
#include "HashPtrTable.H"
#include "fvcDiv.H"
#include "fvmSup.H"
#include "fvMatrix.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::
HeatAndMassTransferPhaseSystem
(
const fvMesh& mesh
)
:
BasePhaseSystem(mesh)
{
this->generatePairsAndSubModels
(
"heatTransfer",
heatTransferModels_
);
this->generatePairsAndSubModels
(
"massTransfer",
massTransferModels_
);
this->generatePairsAndSubModels
(
"interfaceComposition",
interfaceCompositionModels_
);
forAllConstIter
(
phaseSystem::phasePairTable,
this->phasePairs_,
phasePairIter
)
{
const phasePair& pair(phasePairIter());
if (pair.ordered())
{
continue;
}
// Initialy assume no mass transfer
dmdt_.insert
(
pair,
new volScalarField
(
IOobject
(
IOobject::groupName("dmdt", pair.name()),
this->mesh().time().timeName(),
this->mesh(),
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
this->mesh(),
dimensionedScalar("zero", dimDensity/dimTime, 0)
)
);
dmdtExplicit_.insert
(
pair,
new volScalarField
(
IOobject
(
IOobject::groupName("dmdtExplicit", pair.name()),
this->mesh().time().timeName(),
this->mesh()
),
this->mesh(),
dimensionedScalar("zero", dimDensity/dimTime, 0)
)
);
volScalarField H1(heatTransferModels_[pair][pair.first()]->K());
volScalarField H2(heatTransferModels_[pair][pair.second()]->K());
Tf_.insert
(
pair,
new volScalarField
(
IOobject
(
IOobject::groupName("Tf", pair.name()),
this->mesh().time().timeName(),
this->mesh(),
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
(
H1*pair.phase1().thermo().T()
+ H2*pair.phase2().thermo().T()
)
/max
(
H1 + H2,
dimensionedScalar("small", heatTransferModel::dimK, SMALL)
),
zeroGradientFvPatchScalarField::typeName
Henry Weller
committed
)
);
Tf_[pair]->correctBoundaryConditions();
Henry Weller
committed
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::
~HeatAndMassTransferPhaseSystem()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
template<class BasePhaseSystem>
Foam::tmp<Foam::volScalarField>
Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::dmdt
(
const phasePairKey& key
) const
{
const scalar dmdtSign(Pair<word>::compare(dmdt_.find(key).key(), key));
return dmdtSign**dmdt_[key];
}
template<class BasePhaseSystem>
Foam::autoPtr<Foam::phaseSystem::momentumTransferTable>
Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::momentumTransfer() const
Henry Weller
committed
{
autoPtr<phaseSystem::momentumTransferTable>
eqnsPtr(BasePhaseSystem::momentumTransfer());
Henry Weller
committed
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
phaseSystem::momentumTransferTable& eqns = eqnsPtr();
// Source term due to mass trasfer
forAllConstIter
(
phaseSystem::phasePairTable,
this->phasePairs_,
phasePairIter
)
{
const phasePair& pair(phasePairIter());
if (pair.ordered())
{
continue;
}
const volVectorField& U1(pair.phase1().U());
const volVectorField& U2(pair.phase2().U());
const volScalarField dmdt(this->dmdt(pair));
const volScalarField dmdt12(dmdt*pos(dmdt));
const volScalarField dmdt21(dmdt*neg(dmdt));
*eqns[pair.phase1().name()] += fvm::Sp(dmdt21, U1) - dmdt21*U2;
*eqns[pair.phase2().name()] += dmdt12*U1 - fvm::Sp(dmdt12, U2);
}
return eqnsPtr;
}
template<class BasePhaseSystem>
Foam::autoPtr<Foam::phaseSystem::heatTransferTable>
Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::heatTransfer() const
{
autoPtr<phaseSystem::heatTransferTable> eqnsPtr
(
new phaseSystem::heatTransferTable()
);
phaseSystem::heatTransferTable& eqns = eqnsPtr();
forAllConstIter
(
phaseSystem::phaseModelTable,
this->phaseModels_,
phaseModelIter
)
{
const phaseModel& phase(phaseModelIter());
eqns.insert
(
phase.name(),
new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
);
}
// Heat transfer with the interface
forAllConstIter
(
heatTransferModelTable,
heatTransferModels_,
heatTransferModelIter
)
{
const phasePair& pair
(
this->phasePairs_[heatTransferModelIter.key()]
);
const phaseModel* phase = &pair.phase1();
const phaseModel* otherPhase = &pair.phase2();
const volScalarField& Tf(*Tf_[pair]);
const volScalarField K1
(
heatTransferModelIter()[pair.first()]->K()
);
const volScalarField K2
(
heatTransferModelIter()[pair.second()]->K()
);
const volScalarField KEff
(
K1*K2
/max
(
K1 + K2,
dimensionedScalar("small", heatTransferModel::dimK, SMALL)
)
);
const volScalarField* K = &K1;
const volScalarField* otherK = &K2;
forAllConstIter(phasePair, pair, iter)
{
const volScalarField& he(phase->thermo().he());
volScalarField Cpv(phase->thermo().Cpv());
*eqns[phase->name()] +=
(*K)*(Tf - phase->thermo().T())
+ KEff/Cpv*he - fvm::Sp(KEff/Cpv, he);
Swap(phase, otherPhase);
Swap(K, otherK);
}
}
Henry Weller
committed
// Source term due to mass trasfer
forAllConstIter
(
phaseSystem::phasePairTable,
this->phasePairs_,
phasePairIter
)
{
const phasePair& pair(phasePairIter());
if (pair.ordered())
{
continue;
}
const phaseModel& phase1 = pair.phase1();
const phaseModel& phase2 = pair.phase2();
const volScalarField& he1(phase1.thermo().he());
const volScalarField& he2(phase2.thermo().he());
Henry Weller
committed
const volScalarField& K1(phase1.K());
const volScalarField& K2(phase2.K());
Henry Weller
committed
const volScalarField dmdt(this->dmdt(pair));
const volScalarField dmdt12(dmdt*pos(dmdt));
const volScalarField dmdt21(dmdt*neg(dmdt));
const volScalarField& Tf(*Tf_[pair]);
Henry Weller
committed
*eqns[phase1.name()] +=
fvm::Sp(dmdt21, he1) + dmdt21*K1
- dmdt21*(phase2.thermo().he(phase2.thermo().p(), Tf) + K2);
Henry Weller
committed
*eqns[phase2.name()] +=
dmdt12*(phase1.thermo().he(phase1.thermo().p(), Tf) + K1)
- fvm::Sp(dmdt12, he2) - dmdt12*K2;
Henry Weller
committed
}
Henry Weller
committed
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
return eqnsPtr;
}
template<class BasePhaseSystem>
Foam::autoPtr<Foam::phaseSystem::massTransferTable>
Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::massTransfer() const
{
// Create a mass transfer matrix for each species of each phase
autoPtr<phaseSystem::massTransferTable> eqnsPtr
(
new phaseSystem::massTransferTable()
);
phaseSystem::massTransferTable& eqns = eqnsPtr();
forAllConstIter
(
phaseSystem::phaseModelTable,
this->phaseModels_,
phaseModelIter
)
{
const phaseModel& phase(phaseModelIter());
const PtrList<volScalarField>& Yi = phase.Y();
forAll(Yi, i)
{
eqns.insert
(
Yi[i].name(),
new fvScalarMatrix(Yi[i], dimMass/dimTime)
);
}
}
// Reset the interfacial mass flow rates
forAllConstIter
(
phaseSystem::phasePairTable,
this->phasePairs_,
phasePairIter
)
{
const phasePair& pair(phasePairIter());
if (pair.ordered())
{
continue;
}
*dmdt_[pair] =
*dmdtExplicit_[pair];
*dmdtExplicit_[pair] =
dimensionedScalar("zero", dimDensity/dimTime, 0);
}
// Sum up the contribution from each interface composition model
forAllConstIter
(
interfaceCompositionModelTable,
interfaceCompositionModels_,
interfaceCompositionModelIter
)
{
const interfaceCompositionModel& compositionModel
(
interfaceCompositionModelIter()
);
const phasePair& pair
(
this->phasePairs_[interfaceCompositionModelIter.key()]
);
const phaseModel& phase = pair.phase1();
const phaseModel& otherPhase = pair.phase2();
const phasePairKey key(phase.name(), otherPhase.name());
const volScalarField& Tf(*Tf_[key]);
volScalarField& dmdtExplicit(*dmdtExplicit_[key]);
volScalarField& dmdt(*dmdt_[key]);
scalar dmdtSign(Pair<word>::compare(dmdt_.find(key).key(), key));
const volScalarField K
(
massTransferModels_[key][phase.name()]->K()
);
forAllConstIter
(
hashedWordList,
compositionModel.species(),
memberIter
)
{
const word& member = *memberIter;
const word name
(
IOobject::groupName(member, phase.name())
);
const word otherName
(
IOobject::groupName(member, otherPhase.name())
);
const volScalarField KD
(
K*compositionModel.D(member)
);
const volScalarField Yf
(
compositionModel.Yf(member, Tf)
);
// Implicit transport through the phase
*eqns[name] +=
phase.rho()*KD*Yf
- fvm::Sp(phase.rho()*KD, eqns[name]->psi());
// Sum the mass transfer rate
dmdtExplicit += dmdtSign*phase.rho()*KD*Yf;
dmdt -= dmdtSign*phase.rho()*KD*eqns[name]->psi();
// Explicit transport out of the other phase
if (eqns.found(otherName))
{
*eqns[otherName] -=
otherPhase.rho()*KD*compositionModel.dY(member, Tf);
}
}
}
return eqnsPtr;
}
template<class BasePhaseSystem>
void Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::correctThermo()
{
BasePhaseSystem::correctThermo();
// This loop solves for the interface temperatures, Tf, and updates the
// interface composition models.
//
// The rate of heat transfer to the interface must equal the latent heat
// consumed at the interface, i.e.:
//
// H1*(T1 - Tf) + H2*(T2 - Tf) == mDotL
// == K*rho*(Yfi - Yi)*Li
//
// Yfi is likely to be a strong non-linear (typically exponential) function
// of Tf, so the solution for the temperature is newton-accelerated
forAllConstIter
(
phaseSystem::phasePairTable,
this->phasePairs_,
phasePairIter
)
{
const phasePair& pair(phasePairIter());
if (pair.ordered())
{
continue;
}
const phasePairKey key12(pair.first(), pair.second(), true);
const phasePairKey key21(pair.second(), pair.first(), true);
volScalarField H1(heatTransferModels_[pair][pair.first()]->K());
volScalarField H2(heatTransferModels_[pair][pair.second()]->K());
dimensionedScalar HSmall("small", heatTransferModel::dimK, SMALL);
volScalarField mDotL
(
IOobject
(
"mDotL",
this->mesh().time().timeName(),
this->mesh()
),
this->mesh(),
dimensionedScalar("zero", dimEnergy/dimVolume/dimTime, 0)
);
volScalarField mDotLPrime
(
IOobject
(
"mDotLPrime",
this->mesh().time().timeName(),
this->mesh()
),
this->mesh(),
dimensionedScalar("zero", mDotL.dimensions()/dimTemperature, 0)
);
volScalarField& Tf = *Tf_[pair];
// Add latent heats from forward and backward models
if (interfaceCompositionModels_.found(key12))
{
interfaceCompositionModels_[key12]->addMDotL
(
massTransferModels_[pair][pair.first()]->K(),
Tf,
mDotL,
mDotLPrime
);
}
if (interfaceCompositionModels_.found(key21))
{
interfaceCompositionModels_[key21]->addMDotL
(
massTransferModels_[pair][pair.second()]->K(),
Tf,
mDotL,
mDotLPrime
);
}
// Update the interface temperature by applying one step of newton's
// method to the interface relation
Tf -=
(
H1*(Tf - pair.phase1().thermo().T())
+ H2*(Tf - pair.phase2().thermo().T())
+ mDotL
)
/(
max(H1 + H2 + mDotLPrime, HSmall)
);
// Update the interface compositions
if (interfaceCompositionModels_.found(key12))
{
interfaceCompositionModels_[key12]->update(Tf);
}
if (interfaceCompositionModels_.found(key21))
{
interfaceCompositionModels_[key21]->update(Tf);
}
Tf.correctBoundaryConditions();
Henry Weller
committed
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
Info<< "Tf." << pair.name()
<< ": min = " << min(Tf.internalField())
<< ", mean = " << average(Tf.internalField())
<< ", max = " << max(Tf.internalField())
<< endl;
}
}
template<class BasePhaseSystem>
bool Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::read()
{
if (BasePhaseSystem::read())
{
bool readOK = true;
// Models ...
return readOK;
}
else
{
return false;
}
}
// ************************************************************************* //