alphatWallBoilingWallFunctionFvPatchScalarField.C 43.4 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
5
    \\  /    A nd           | www.openfoam.com
6
7
     \\/     M anipulation  |
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
8
9
    Copyright (C) 2015-2018 OpenFOAM Foundation
    Copyright (C) 2018 OpenCFD Ltd
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
-------------------------------------------------------------------------------
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 "alphatWallBoilingWallFunctionFvPatchScalarField.H"
#include "fvPatchFieldMapper.H"
#include "addToRunTimeSelectionTable.H"

#include "phaseSystem.H"
#include "compressibleTurbulenceModel.H"
#include "ThermalDiffusivity.H"
#include "PhaseCompressibleTurbulenceModel.H"
#include "saturationModel.H"
#include "wallFvPatch.H"
#include "uniformDimensionedFields.H"
#include "mathematicalConstants.H"

using namespace Foam::constant::mathematical;

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

const Foam::Enum
<
    Foam::compressible::
    alphatWallBoilingWallFunctionFvPatchScalarField::phaseType
>
Foam::compressible::
alphatWallBoilingWallFunctionFvPatchScalarField::phaseTypeNames_
{
    { phaseType::vaporPhase, "vapor" },
    { phaseType::liquidPhase, "liquid" },
};


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{
namespace compressible
{

// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

alphatWallBoilingWallFunctionFvPatchScalarField::
alphatWallBoilingWallFunctionFvPatchScalarField
(
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF
)
:
    alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(p, iF),
    otherPhaseName_("vapor"),
    phaseType_(liquidPhase),
78
    relax_(),
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
    AbyV_(p.size(), 0),
    alphatConv_(p.size(), 0),
    dDep_(p.size(), 1e-5),
    qq_(p.size(), 0),
    K_(4),
    partitioningModel_(nullptr),
    nucleationSiteModel_(nullptr),
    departureDiamModel_(nullptr),
    departureFreqModel_(nullptr),
    filmBoilingModel_(nullptr),
    LeidenfrostModel_(nullptr),
    CHFModel_(nullptr),
    CHFSoobModel_(nullptr),
    MHFModel_(nullptr),
    TDNBModel_(nullptr),
    wp_(1)
{
    AbyV_ = this->patch().magSf();
    forAll(AbyV_, facei)
    {
        const label faceCelli = this->patch().faceCells()[facei];
        AbyV_[facei] /= iF.mesh().V()[faceCelli];
    }
}


alphatWallBoilingWallFunctionFvPatchScalarField::
alphatWallBoilingWallFunctionFvPatchScalarField
(
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF,
    const dictionary& dict
)
:
    alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(p, iF, dict),
    otherPhaseName_(dict.lookup("otherPhase")),
    phaseType_(phaseTypeNames_.read(dict.lookup("phaseType"))),
116
    relax_(Function1<scalar>::New("relax", dict)),
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
142
143
144
145
146
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
176
177
178
179
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
    AbyV_(p.size(), 0),
    alphatConv_(p.size(), 0),
    dDep_(p.size(), 1e-5),
    qq_(p.size(), 0),
    K_(4),
    partitioningModel_(nullptr),
    nucleationSiteModel_(nullptr),
    departureDiamModel_(nullptr),
    departureFreqModel_(nullptr),
    filmBoilingModel_(nullptr),
    LeidenfrostModel_(nullptr),
    CHFModel_(nullptr),
    CHFSoobModel_(nullptr),
    MHFModel_(nullptr),
    TDNBModel_(nullptr),
    wp_(1)
{

    // Check that otherPhaseName != this phase
    if (internalField().group() == otherPhaseName_)
    {
        FatalErrorInFunction
            << "otherPhase should be the name of the vapor phase that "
            << "corresponds to the liquid base of vice versa" << nl
            << "This phase: " << internalField().group() << nl
            << "otherPhase: " << otherPhaseName_
            << abort(FatalError);
    }

    switch (phaseType_)
    {
        case vaporPhase:
        {
            partitioningModel_ =
                wallBoilingModels::partitioningModel::New
                (
                    dict.subDict("partitioningModel")
                );

            dmdt_ = 0;

            break;
        }
        case liquidPhase:
        {
            partitioningModel_ =
                wallBoilingModels::partitioningModel::New
                (
                    dict.subDict("partitioningModel")
                );

            nucleationSiteModel_ =
                wallBoilingModels::nucleationSiteModel::New
                (
                    dict.subDict("nucleationSiteModel")
                );

            departureDiamModel_ =
                wallBoilingModels::departureDiameterModel::New
                (
                    dict.subDict("departureDiamModel")
                );

            departureFreqModel_ =
                wallBoilingModels::departureFrequencyModel::New
                (
                    dict.subDict("departureFreqModel")
                );

            const dictionary* LeidenfrostDict =
                dict.findDict("LeidenfrostModel");

            if (LeidenfrostDict)
            {
                LeidenfrostModel_ =
                    wallBoilingModels::LeidenfrostModel::New(*LeidenfrostDict);
            }

            const dictionary* CHFDict = dict.findDict("CHFModel");

            if (CHFDict)
            {
                CHFModel_ =
                    wallBoilingModels::CHFModel::New(*CHFDict);
            }

            const dictionary* HFSubCoolDict = dict.findDict("CHFSubCoolModel");

            if (HFSubCoolDict)
            {
                CHFSoobModel_ =
                    wallBoilingModels::CHFSubCoolModel::New(*HFSubCoolDict);
            }

            const dictionary* MHFDict = dict.findDict("MHFModel");

            if (MHFDict)
            {
                MHFModel_ =
                    wallBoilingModels::MHFModel::New(*MHFDict);
            }

            const dictionary* TDNBDict = dict.findDict("TDNBModel");

            if (TDNBDict)
            {
                TDNBModel_ =
                    wallBoilingModels::TDNBModel::New(*TDNBDict);
            }

            const dictionary* filmDict = dict.findDict("filmBoilingModel");

            if (filmDict)
            {
                filmBoilingModel_ =
                    wallBoilingModels::filmBoilingModel::New(*filmDict);
            }

            if (dict.found("dDep"))
            {
                dDep_ = scalarField("dDep", dict, p.size());
            }

            if (dict.found("K"))
            {
                dict.lookup("K") >> K_;
            }

            if (dict.found("wp"))
            {
                dict.lookup("wp") >> wp_;
            }

            if (dict.found("qQuenching"))
            {
                qq_ = scalarField("qQuenching", dict, p.size());
            }

            break;
        }
    }

    if (dict.found("alphatConv"))
    {
        alphatConv_ = scalarField("alphatConv", dict, p.size());
    }

    AbyV_ = this->patch().magSf();
    forAll(AbyV_, facei)
    {
        const label faceCelli = this->patch().faceCells()[facei];
        AbyV_[facei] /= iF.mesh().V()[faceCelli];
    }
}


alphatWallBoilingWallFunctionFvPatchScalarField::
alphatWallBoilingWallFunctionFvPatchScalarField
(
    const alphatWallBoilingWallFunctionFvPatchScalarField& psf,
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF,
    const fvPatchFieldMapper& mapper
)
:
    alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
    (
        psf,
        p,
        iF,
        mapper
    ),
    otherPhaseName_(psf.otherPhaseName_),
    phaseType_(psf.phaseType_),
291
    relax_(psf.relax_.clone()),
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
    AbyV_(psf.AbyV_),
    alphatConv_(psf.alphatConv_, mapper),
    dDep_(psf.dDep_, mapper),
    qq_(psf.qq_, mapper),
    K_(psf.K_),
    partitioningModel_(psf.partitioningModel_),
    nucleationSiteModel_(psf.nucleationSiteModel_),
    departureDiamModel_(psf.departureDiamModel_),
    filmBoilingModel_(psf.filmBoilingModel_),
    LeidenfrostModel_(psf.LeidenfrostModel_),
    CHFModel_(psf.CHFModel_),
    CHFSoobModel_(psf.CHFSoobModel_),
    MHFModel_(psf.MHFModel_),
    TDNBModel_(psf.TDNBModel_),
    wp_(psf.wp_)
{}


alphatWallBoilingWallFunctionFvPatchScalarField::
alphatWallBoilingWallFunctionFvPatchScalarField
(
    const alphatWallBoilingWallFunctionFvPatchScalarField& psf
)
:
    alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(psf),
    otherPhaseName_(psf.otherPhaseName_),
    phaseType_(psf.phaseType_),
319
    relax_(psf.relax_.clone()),
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
    AbyV_(psf.AbyV_),
    alphatConv_(psf.alphatConv_),
    dDep_(psf.dDep_),
    qq_(psf.qq_),
    K_(psf.K_),
    partitioningModel_(psf.partitioningModel_),
    nucleationSiteModel_(psf.nucleationSiteModel_),
    departureDiamModel_(psf.departureDiamModel_),
    filmBoilingModel_(psf.filmBoilingModel_),
    LeidenfrostModel_(psf.LeidenfrostModel_),
    CHFModel_(psf.CHFModel_),
    CHFSoobModel_(psf.CHFSoobModel_),
    MHFModel_(psf.MHFModel_),
    TDNBModel_(psf.TDNBModel_),
    wp_(psf.wp_)
{}


alphatWallBoilingWallFunctionFvPatchScalarField::
alphatWallBoilingWallFunctionFvPatchScalarField
(
    const alphatWallBoilingWallFunctionFvPatchScalarField& psf,
    const DimensionedField<scalar, volMesh>& iF
)
:
    alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(psf, iF),
    otherPhaseName_(psf.otherPhaseName_),
    phaseType_(psf.phaseType_),
348
    relax_(psf.relax_.clone()),
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
    AbyV_(psf.AbyV_),
    alphatConv_(psf.alphatConv_),
    dDep_(psf.dDep_),
    qq_(psf.qq_),
    K_(psf.K_),
    partitioningModel_(psf.partitioningModel_),
    nucleationSiteModel_(psf.nucleationSiteModel_),
    departureDiamModel_(psf.departureDiamModel_),
    filmBoilingModel_(psf.filmBoilingModel_),
    LeidenfrostModel_(psf.LeidenfrostModel_),
    CHFModel_(psf.CHFModel_),
    CHFSoobModel_(psf.CHFSoobModel_),
    MHFModel_(psf.MHFModel_),
    TDNBModel_(psf.TDNBModel_),
    wp_(psf.wp_)
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

bool alphatWallBoilingWallFunctionFvPatchScalarField::
activePhasePair(const phasePairKey& phasePair) const
{
    if (phasePair == phasePairKey(otherPhaseName_, internalField().group()))
    {
        return true;
    }
    else
    {
        return false;
    }
}

const scalarField& alphatWallBoilingWallFunctionFvPatchScalarField::
dmdt(const phasePairKey& phasePair) const
{
    if (activePhasePair(phasePair))
    {
        return dmdt_;
    }
    else
    {
        FatalErrorInFunction
            << " dmdt requested for invalid phasePair!"
            << abort(FatalError);

        return dmdt_;
    }
}

const scalarField& alphatWallBoilingWallFunctionFvPatchScalarField::
mDotL(const phasePairKey& phasePair) const
{
    if (activePhasePair(phasePair))
    {
        return mDotL_;
    }
    else
    {
        FatalErrorInFunction
            << " mDotL requested for invalid phasePair!"
            << abort(FatalError);

        return mDotL_;
    }
}

void alphatWallBoilingWallFunctionFvPatchScalarField::updateCoeffs()
{

    if (updated())
    {
        return;
    }

    // Check that partitioningModel has been constructed
    if (!partitioningModel_.valid())
    {
        FatalErrorInFunction
            << "partitioningModel has not been constructed!"
            << abort(FatalError);
    }

    // Lookup the fluid model
    const phaseSystem& fluid =
        refCast<const phaseSystem>
        (
            db().lookupObject<phaseSystem>("phaseProperties")
        );

    const saturationModel& satModel =
        db().lookupObject<saturationModel>("saturationModel");

    const label patchi = patch().index();

444
445
446
    const scalar t = this->db().time().timeOutputValue();
    scalar relax = relax_->value(t);

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
    switch (phaseType_)
    {
        case vaporPhase:
        {
            const phaseModel& vapor
            (
                fluid.phases()[internalField().group()]
            );

            const fvPatchScalarField& hewv =
                vapor.thermo().he().boundaryField()[patchi];

            // Vapor Liquid phase fraction at the wall
            const scalarField vaporw(vapor.boundaryField()[patchi]);

            // NOTE! Assumes 1-thisPhase for liquid fraction in
            // multiphase simulations
            const scalarField fLiquid
            (
                partitioningModel_->fLiquid(1-vaporw)
            );

            const tmp<scalarField> talphaw = vapor.thermo().alpha(patchi);
            const scalarField& alphaw = talphaw();

            const scalarField heSnGrad(max(hewv.snGrad(), scalar(1e-16)));

            // Convective thermal diffusivity for single phase
            const scalarField alphatv(calcAlphat(*this));

            forAll (*this, i)
            {
479
480
481
482
483
                this->operator[](i) =
                (
                    (1 - fLiquid[i])*(alphatv[i])
                   /max(vaporw[i], scalar(1e-8))
                );
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
581
582
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
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
            }

            if (debug)
            {
                Info<< "alphat for vapour : " << nl << endl;

                Info<< "  alphatEffv: " << gMin(vaporw*(*this + alphaw))
                    << " - " << gMax(vaporw*(*this + alphaw)) << endl;

                const scalarField qEff(vaporw*(*this + alphaw)*hewv.snGrad());

                scalar Qeff = gSum(qEff*patch().magSf());
                Info<< " Effective heat transfer rate to vapor:" << Qeff
                    << nl << endl;
            }
            break;
        }
        case liquidPhase:
        {
            // Check that nucleationSiteModel has been constructed
            if (!nucleationSiteModel_.valid())
            {
                FatalErrorInFunction
                    << "nucleationSiteModel has not been constructed!"
                    << abort(FatalError);
            }

            // Check that departureDiameterModel has been constructed
            if (!departureDiamModel_.valid())
            {
                FatalErrorInFunction
                    << "departureDiameterModel has not been constructed!"
                    << abort(FatalError);
            }

            // Check that nucleationSiteModel has been constructed
            if (!departureFreqModel_.valid())
            {
                FatalErrorInFunction
                    << "departureFrequencyModel has not been constructed!"
                    << abort(FatalError);
            }

            const phaseModel& liquid
            (
                fluid.phases()[internalField().group()]
            );

            const phaseModel& vapor(fluid.phases()[otherPhaseName_]);

            // Retrieve turbulence properties from models
            const phaseCompressibleTurbulenceModel& turbModel =
                db().lookupObject<phaseCompressibleTurbulenceModel>
                (
                    IOobject::groupName
                    (
                        turbulenceModel::propertiesName,
                        liquid.name()
                    )
                );
            const phaseCompressibleTurbulenceModel& vaporTurbModel =
                db().lookupObject<phaseCompressibleTurbulenceModel>
                (
                    IOobject::groupName
                    (
                        turbulenceModel::propertiesName,
                        vapor.name()
                    )
                );

            const tmp<scalarField> tnutw = turbModel.nut(patchi);

            const scalar Cmu25(pow025(Cmu_));

            const scalarField& y = turbModel.y()[patchi];

            const tmp<scalarField> tmuw = turbModel.mu(patchi);
            const scalarField& muw = tmuw();

            const tmp<scalarField> talphaw = liquid.thermo().alphahe(patchi);
            const scalarField& alphaw = talphaw();

            const tmp<volScalarField> tk = turbModel.k();
            const volScalarField& k = tk();
            const fvPatchScalarField& kw = k.boundaryField()[patchi];

            const fvPatchVectorField& Uw =
                turbModel.U().boundaryField()[patchi];
            const scalarField magUp(mag(Uw.patchInternalField() - Uw));
            const scalarField magGradUw(mag(Uw.snGrad()));

            const fvPatchScalarField& rhow =
                turbModel.rho().boundaryField()[patchi];


            const fvPatchScalarField& Tw =
                liquid.thermo().T().boundaryField()[patchi];
            const scalarField Tc(Tw.patchInternalField());

            const scalarField uTau(Cmu25*sqrt(kw));

            const scalarField yPlus(uTau*y/(muw/rhow));

            const scalarField Pr(muw/alphaw);

            // Molecular-to-turbulent Prandtl number ratio
            const scalarField Prat(Pr/Prt_);

            // Thermal sublayer thickness
            const scalarField P(this->Psmooth(Prat));

            const scalarField yPlusTherm(this->yPlusTherm(P, Prat));

            const fvPatchScalarField& rhoVaporw =
                vaporTurbModel.rho().boundaryField()[patchi];

            tmp<volScalarField> tCp = liquid.thermo().Cp();
            const volScalarField& Cp = tCp();
            const fvPatchScalarField& Cpw = Cp.boundaryField()[patchi];

            // Saturation temperature
            const tmp<volScalarField> tTsat =
                satModel.Tsat(liquid.thermo().p());

            const volScalarField& Tsat = tTsat();
            const fvPatchScalarField& Tsatw(Tsat.boundaryField()[patchi]);
            const scalarField Tsatc(Tsatw.patchInternalField());

            const fvPatchScalarField& pw =
                liquid.thermo().p().boundaryField()[patchi];

            const fvPatchScalarField& hew =
                liquid.thermo().he().boundaryField()[patchi];

            const scalarField hw
            (
                liquid.thermo().he().member() == "e"
              ? hew.patchInternalField() + pw/rhow.patchInternalField()
              : hew.patchInternalField()
            );

            const scalarField L
            (
                vapor.thermo().he().member() == "e"
              ? vapor.thermo().he(pw, Tsatc, patchi) + pw/rhoVaporw - hw
              : vapor.thermo().he(pw, Tsatc, patchi) - hw
            );

            // Liquid phase fraction at the wall
            const scalarField liquidw(liquid.boundaryField()[patchi]);

            const scalarField fLiquid(partitioningModel_->fLiquid(liquidw));

637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
            // Liquid temperature at y+=250 is estimated from logarithmic
            // thermal wall function (Koncar, Krepper & Egorov, 2005)
            const scalarField Tplus_y250
            (
                Prt_*(log(E_*250)/kappa_ + P)
            );
            const scalarField Tplus(Prt_*(log(E_*yPlus)/kappa_ + P));
            scalarField Tl(Tw - (Tplus_y250/Tplus)*(Tw - Tc));
            Tl = max(Tc - 40, Tl);

            // Film, transient boiling regimes
            scalarField Qtb(this->size(), 0);
            scalarField tDNB(this->size(), GREAT);
            scalarField TLeiden(this->size(), GREAT);
            scalarField htcFilmBoiling(this->size(), 0);

            if
            (
                CHFModel_.valid()
                && CHFSoobModel_.valid()
                && TDNBModel_.valid()
                && MHFModel_.valid()
                && LeidenfrostModel_.valid()
                && filmBoilingModel_.valid()
            )
662
            {
663
664

                const scalarField CHF
665
                (
666
667
668
669
670
671
672
673
674
675
                    CHFModel_->CHF
                    (
                        liquid,
                        vapor,
                        patchi,
                        Tl,
                        Tsatw,
                        L
                    )
                );
676

677
678
679
680
                // Effect of sub-cooling to the CHF in saturated conditions
                const scalarField CHFSubCool
                (
                    CHFSoobModel_->CHFSubCool
681
                    (
682
683
684
685
686
687
688
689
690
691
                        liquid,
                        vapor,
                        patchi,
                        Tl,
                        Tsatw,
                        L
                    )
                );

                const scalarField CHFtotal(CHF*CHFSubCool);
692

693
694
                tDNB =
                    TDNBModel_->TDNB
695
                    (
696
697
698
699
700
701
                        liquid,
                        vapor,
                        patchi,
                        Tl,
                        Tsatw,
                        L
702
703
                    );

704
705
706
707
708
709
710
711
712
713
714
715
                const scalarField MHF
                (
                    MHFModel_->MHF
                    (
                        liquid,
                        vapor,
                        patchi,
                        Tl,
                        Tsatw,
                        L
                    )
                );
716

717
718
719
720
721
722
723
724
725
726
                TLeiden =
                    LeidenfrostModel_->TLeid
                    (
                        liquid,
                        vapor,
                        patchi,
                        Tl,
                        Tsatw,
                        L
                    );
727

728
729
730
                // htc for film boiling
                htcFilmBoiling =
                    filmBoilingModel_->htcFilmBoil
731
                    (
732
733
734
735
736
737
                        liquid,
                        vapor,
                        patchi,
                        Tl,
                        Tsatw,
                        L
738
739
                    );

740
741
742
743
744
745
746
747
                // htc for film transition boiling

                // Indicator between CHF (phi = 0) and MHF (phi = 1)
                const scalarField phi
                (
                    min
                    (
                        max
748
                        (
749
750
751
752
753
754
                            wp_*(Tw - tDNB)/(TLeiden - tDNB),
                            scalar(0)
                        ),
                        scalar(1)
                    )
                );
755

756
757
758
759
760
761
762
763
                Qtb = CHFtotal*(1 - phi) + phi*MHF;


                if (debug & 2)
                {

                    scalarField& phiFluid =
                        const_cast<fvPatchScalarField&>
764
                        (
765
766
767
768
769
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("phiTb")
770
771
                        );

772
                    phiFluid = phi;
773

774
775
                    scalarField& CHFFluid =
                        const_cast<fvPatchScalarField&>
776
                        (
777
778
779
780
781
782
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("CHFtotal")
                        );
783

784
                    CHFFluid = CHFtotal;
785
786
                }

787
            }
788
789


790
791
792
793
            // Sub-cool boiling Nucleation
            const scalarField N
            (
                nucleationSiteModel_->N
794
795
796
797
798
799
800
                (
                    liquid,
                    vapor,
                    patchi,
                    Tl,
                    Tsatw,
                    L
801
802
803
804
805
806
807
808
809
810
811
812
813
                )
            );

            // Bubble departure diameter:
            dDep_ = departureDiamModel_->dDeparture
            (
                liquid,
                vapor,
                patchi,
                Tl,
                Tsatw,
                L
            );
814

815
816
817
818
            // Bubble departure frequency:
            const scalarField fDep
            (
                departureFreqModel_->fDeparture
819
                (
820
821
822
823
824
825
                    liquid,
                    vapor,
                    patchi,
                    dDep_
                )
            );
826

827
828
            // Convective thermal diffusivity for single phase
            alphatConv_ = calcAlphat(alphatConv_);
829

830
831
            // Convective heat transfer area for Sub-cool boiling
            scalarField A1(this->size(), 0);
832

833
            const scalarField hewSn(hew.snGrad());
834

835
            scalarField alphaFilm(this->size(), 0);
836

837
838
            // Use to identify regimes per face
            labelField regimeTypes(A1.size(), -1);
839

840
841
842
            forAll (*this, i)
            {
                if (Tw[i] > Tsatw[i])
843
                {
844
845
                    // Sub-cool boiling
                    if (Tw[i] < tDNB[i])
846
847
                    {
                        // Sub-cool boiling
848
                        regimeTypes[i] = regimeType::subcool;
849

850
851
852
853
                        // Del Valle & Kenning (1985)
                        const scalar Ja =
                            rhow[i]*Cpw[i]*(Tsatw[i] - Tl[i])
                            /(rhoVaporw[i]*L[i]);
854

855
856
                        const scalar Al =
                            fLiquid[i]*4.8*exp(min(-Ja/80, log(VGREAT)));
857

858
859
860
861
                        const scalar A2
                        (
                            min(pi*sqr(dDep_[i])*N[i]*Al/4, scalar(1))
                        );
862

863
                        A1[i] = max(1 - A2, scalar(1e-4));
864

865
866
867
868
                        // Following Bowring(1962)
                        const scalar A2E
                        (
                            min
869
                            (
870
871
872
873
                                pi*sqr(dDep_[i])*N[i]*Al/4,
                                scalar(5)
                            )
                        );
874

875
876
877
                        // Volumetric mass source in the near wall cell due
                        // to the wall boiling
                        dmdt_[i] =
878
                            (
879
880
881
                                (1 - relax)*dmdt_[i]
                                + relax*(1.0/6.0)*A2E*dDep_[i]*rhoVaporw[i]
                                * fDep[i]*AbyV_[i]
882
883
                            );

884
885
886
                        // Volumetric source in the near wall cell due to
                        // the wall boiling
                        mDotL_[i] = dmdt_[i]*L[i];
887

888
889
890
891
892
                        // Quenching heat transfer coefficient
                        const scalar hQ
                        (
                            2*(alphaw[i]*Cpw[i])*fDep[i]
                            *sqrt
893
                            (
894
895
896
                                (0.8/max(fDep[i], SMALL))/(pi*alphaw[i]/rhow[i])
                            )
                        );
897

898
899
                        // Quenching heat flux in Sub-cool boiling
                        qq_[i] =
900
                            (
901
902
                                (1 - relax)*qq_[i]
                                + relax*A2*hQ*max(Tw[i] - Tl[i], scalar(0))
903
904
                            );

905
906
907
                        this->operator[](i) =
                        (
                            max
908
                            (
909
910
911
912
                                A1[i]*alphatConv_[i]
                                + (
                                    (qq_[i] + mDotL_[i]/AbyV_[i])
                                    / max(hewSn[i], scalar(1e-16))
913
                                )
914
915
916
917
918
919
920
921
922
                                /max(liquidw[i], scalar(1e-8)),
                                1e-8
                            )
                        );
                    }
                    else if (Tw[i] > tDNB[i] && Tw[i] < TLeiden[i])
                    {
                        // transient boiling
                        regimeTypes[i] = regimeType::transient;
923

924
925
                        // No convective heat tranfer
                        alphatConv_[i] = 0.0;
926

927
928
                        // transient boiling
                        dmdt_[i] =
929
                                fLiquid[i]
930
931
932
                                *(
                                    relax*Qtb[i]*AbyV_[i]/L[i]
                                    + (1 - relax)*dmdt_[i]
933
934
                                );

935
                        mDotL_[i] = dmdt_[i]*L[i];
936
937


938
939
                        // No quenching flux
                        qq_[i] = 0.0;
940

941
942
943
                        this->operator[](i) =
                        max
                        (
944
                            (
945
946
947
948
949
                                mDotL_[i]/AbyV_[i]
                                /max(hewSn[i], scalar(1e-16))
                            )/max(liquidw[i], scalar(1e-8)),
                            1e-8
                        );
950
951

                    }
952
                    else if (Tw[i] > TLeiden[i])
953
                    {
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
                        regimeTypes[i] = regimeType::film; // film boiling

                        // No convective heat tranfer
                        alphatConv_[i] = 0.0;

                        // Film boiling
                        dmdt_[i] =
                            fLiquid[i]
                            *(
                                relax*htcFilmBoiling[i]
                                *max(Tw[i] - Tsatw[i], 0)
                                *AbyV_[i]/L[i]
                                + (1 - relax)*dmdt_[i]
                            );


                        mDotL_[i] = dmdt_[i]*L[i];

                        // No quenching flux
973
974
                        qq_[i] = 0.0;

975
                        alphaFilm[i] =
976
                        (
977
                            mDotL_[i]/AbyV_[i]/max(hewSn[i], scalar(1e-16))
978
                        );
979
980
981
982
983
984
985
986
987

                        // alphat is added alphal and multiplied by phase
                        // alphaFilm in the coupled BC. We substract
                        // alpha and divide by phase to get a net alphaFilm
                        this->operator[](i) =
                            (
                                alphaFilm[i]/max(liquidw[i], scalar(1e-8))
                              - alphaw[i]
                            );
988
989
                    }
                }
990
991
992
993
994
                else
                {
                    // Tw below Tsat. No boiling single phase convection
                    // Single phase
                    regimeTypes[i] = regimeType::nonBoiling;
995

996
997
998
                    qq_[i] = 0.0;
                    mDotL_[i] = 0.0;
                    dmdt_[i] = 0.0;
999

1000
1001
                    // Turbulente thermal diffusivity for single phase.
                    this->operator[](i) =
1002
                    (
1003
1004
1005
1006
1007
1008
                        max
                        (
                            fLiquid[i]*(alphatConv_[i])
                            /max(liquidw[i], scalar(1e-8)),
                            1e-8
                        )
1009
                    );
1010
1011
                }
            }
1012

1013
1014
1015
1016
1017
1018
            if (debug)
            {
                const scalarField qEff
                (
                    liquidw*(*this + alphaw)*hew.snGrad()
                );
1019

1020
                Info<< "alphat for liquid:  " <<  nl << endl;
1021

1022
1023
                Info<< "  alphatl: " << gMin((*this)) << " - "
                    << gMax((*this)) << endl;
1024

1025
1026
                Info<< "  dmdt: " << gMin((dmdt_)) << " - "
                    << gMax((dmdt_)) << endl;
1027

1028
1029
                Info<< "  alphatlEff: " << gMin(liquidw*(*this + alphaw))
                    << " - " << gMax(liquidw*(*this + alphaw)) << endl;
1030

1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
                scalar Qeff = gSum(qEff*patch().magSf());
                Info<< " Effective heat transfer rate to liquid: " << Qeff
                    << endl << nl;

                if (debug & 2)
                {
                    /*
                    scalarField& faceRegimes =
                        const_cast<fvPatchScalarField&>
                        (
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("faceRegimes")
                        );
                    */

                    scalar nSubCool(0);
                    scalar nTransient(0);
                    scalar nFilm(0);
                    scalar nNonBoiling(0);
1053

1054
1055
1056
1057
                    scalarField nSubCools(this->size(), 0);
                    scalarField nTransients(this->size(), 0);
                    scalarField nFilms(this->size(), 0);
                    scalarField nNonBoilings(this->size(), 0);
1058

1059
1060
1061
1062
                    forAll (*this, i)
                    {
                        //faceRegimes[i] = regimeTypes[i];
                        switch (regimeTypes[i])
1063
                        {
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
                            case regimeType::subcool:
                                nSubCool++;
                                nSubCools[i] = 1;
                            break;

                            case regimeType::transient:
                                nTransient++;
                                nTransients[i] = 1;
                            break;

                            case regimeType::film:
                                nFilm++;
                                nFilms[i] = 1;
                            break;

                            case regimeType::nonBoiling:
                                nNonBoiling++;
                                nNonBoilings[i] = 1;
                            break;
1083
                        }
1084
                    }
1085

1086
                    Info<< "Faces regime :  " <<  nl << endl;
1087

1088
1089
1090
1091
1092
                    Info<< "    sub Cool faces : " << nSubCool << endl;
                    Info<< "    transient faces : " << nTransient << endl;
                    Info<< "    film faces : " << nFilm << endl;
                    Info<< "    non-Boiling faces : " << nNonBoiling << endl;
                    Info<< "    total faces : " << this->size() << endl << nl;
1093

1094
1095
1096
1097
1098
1099
1100
1101
                    const scalarField qc
                    (
                        nNonBoilings*fLiquid*A1*(alphatConv_ + alphaw)
                        *hew.snGrad()
                    );
                    /*
                    scalarField& qcFluid =
                        const_cast<fvPatchScalarField&>
1102
                        (
1103
1104
1105
1106
1107
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("qc")
1108
1109
                        );

1110
1111
1112
1113
1114
                    qcFluid = qc;
                    */

                    scalar Qc = gSum(qc*patch().magSf());
                    Info<< " Convective heat transfer: " << Qc << endl;
1115

1116
1117
1118
1119
1120
1121
1122
                    const scalarField qFilm
                    (
                        relax*fLiquid*nFilms*htcFilmBoiling*(Tw - Tsatw)
                    );
                    /*
                    scalarField& qFilmFluid =
                        const_cast<fvPatchScalarField&>
1123
                        (
1124
1125
1126
1127
1128
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("qFilm")
1129
1130
                        );

1131
1132
                    qFilmFluid = qFilm;
                    */
1133

1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
                    /*
                    scalarField& htcFilm =
                        const_cast<fvPatchScalarField&>
                    (
                        patch().lookupPatchField
                        <
                            volScalarField,
                            scalar
                        >("htcFilmBoiling")
                    );
1144

1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
                    htcFilm = htcFilmBoiling;
                    */
                    scalar QFilm = gSum(qFilm*patch().magSf());
                    Info<< " Film boiling heat transfer: " << QFilm << endl;

                    Info<< " Htc Film Boiling coeff: "
                        << gMin(nFilms*htcFilmBoiling)
                        << " - "
                        << gMax(nFilms*htcFilmBoiling) << endl;
                    /*
                    scalarField& qtbFluid =
                        const_cast<fvPatchScalarField&>
1157
                        (
1158
1159
1160
1161
1162
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("qtb")
1163
                        );
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
                    qtbFluid = fLiquid*Qtb;
                    */
                    
                    scalar Qtbtot =
                        gSum(fLiquid*nTransients*Qtb*patch().magSf());
                    Info<< " Transient boiling heat transfer:" << Qtbtot
                        << endl;

                    /*
                    scalarField& TdnbFluid =
                        const_cast<fvPatchScalarField&>
                        (
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("Tdnb")
                        );
                    TdnbFluid = tDNB;
                    */
1184

1185
1186
                    Info<< " TDNB: " << gMin(tDNB) << " - " << gMax(tDNB)
                        << endl;
1187

1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
                    const scalarField qSubCool
                    (
                        fLiquid*nSubCools*
                        (
                            A1*alphatConv_*hew.snGrad()
                            + qe() + qq()
                        )
                    );
                    
                    /*
                    scalarField& qSubCoolFluid =
                        const_cast<fvPatchScalarField&>
                        (
                            patch().lookupPatchField
                            <
                                volScalarField,
                                scalar
                            >("qSubCool")
                        );
1207

1208
1209
1210
1211
                    qSubCoolFluid = qSubCool;
                    */
                    
                    scalar QsubCool = gSum(qSubCool*patch().magSf());
1212

1213
1214
                    Info<< " Sub Cool boiling heat transfer: " << QsubCool
                        << endl;
1215

1216
1217
1218
1219
1220
1221
1222
1223
                    Info<< "  N: " << gMin(nSubCools*N) << " - "
                        << gMax(nSubCools*N) << endl;
                    Info<< "  dDep: " << gMin(nSubCools*dDep_) << " - "
                        << gMax(nSubCools*dDep_) << endl;
                    Info<< "  fDep: " << gMin(nSubCools*fDep) << " - "
                        << gMax(nSubCools*fDep) << endl;
                    Info<< "  A1: " << gMin(nSubCools*A1) << " - "
                        << gMax(nSubCools*A1) << endl;
1224

1225
                    Info<< nl;
1226
1227
1228
                }
            }
        }
1229
        break;
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
        default:
        {
            FatalErrorInFunction
                << "Unknown phase type. Valid types are: "
                << phaseTypeNames_ << nl << exit(FatalError);
        }
    }

    fixedValueFvPatchScalarField::updateCoeffs();
}


void alphatWallBoilingWallFunctionFvPatchScalarField::write(Ostream& os) const
{
    fvPatchField<scalar>::write(os);

    os.writeKeyword("phaseType") << phaseTypeNames_[phaseType_]
        << token::END_STATEMENT << nl;

1249
    relax_->writeData(os);
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378

    switch (phaseType_)
    {
        case vaporPhase:
        {
            os.writeKeyword("partitioningModel") << nl;
            os  << indent << token::BEGIN_BLOCK << incrIndent << nl;
            partitioningModel_->write(os);
            os << decrIndent << indent << token::END_BLOCK << nl;

            if (filmBoilingModel_.valid())
            {
                os.writeKeyword("filmBoilingModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                filmBoilingModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            if (LeidenfrostModel_.valid())
            {
                os.writeKeyword("LeidenfrostModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                LeidenfrostModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            break;
        }
        case liquidPhase:
        {
            os.writeKeyword("partitioningModel") << nl;
            os << indent << token::BEGIN_BLOCK << incrIndent << nl;
            partitioningModel_->write(os);
            os << decrIndent << indent << token::END_BLOCK << nl;

            os.writeKeyword("nucleationSiteModel") << nl;
            os << indent << token::BEGIN_BLOCK << incrIndent << nl;
            nucleationSiteModel_->write(os);
            os << decrIndent << indent << token::END_BLOCK << nl;

            os.writeKeyword("departureDiamModel") << nl;
            os << indent << token::BEGIN_BLOCK << incrIndent << nl;
            departureDiamModel_->write(os);
            os << decrIndent << indent << token::END_BLOCK << nl;

            os.writeKeyword("departureFreqModel") << nl;
            os << indent << token::BEGIN_BLOCK << incrIndent << nl;
            departureFreqModel_->write(os);
            os << decrIndent << indent << token::END_BLOCK << nl;

            if (filmBoilingModel_.valid())
            {
                os.writeKeyword("filmBoilingModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                filmBoilingModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            if (LeidenfrostModel_.valid())
            {
                os.writeKeyword("LeidenfrostModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                LeidenfrostModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            if (CHFModel_.valid())
            {
                os.writeKeyword("CHFModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                CHFModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            if (CHFSoobModel_.valid())
            {
                os.writeKeyword("CHFSubCoolModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                CHFSoobModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            if (MHFModel_.valid())
            {
                os.writeKeyword("MHFModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                MHFModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            if (TDNBModel_.valid())
            {
                os.writeKeyword("TDNBModel") << nl;
                os << indent << token::BEGIN_BLOCK << incrIndent << nl;
                TDNBModel_->write(os);
                os << decrIndent << indent << token::END_BLOCK << nl;
            }

            os.writeKeyword("K") << K_ << token::END_STATEMENT << nl;
            os.writeKeyword("wp") << wp_ << token::END_STATEMENT << nl;
            break;
        }
    }

    os.writeKeyword("otherPhase") << otherPhaseName_
        << token::END_STATEMENT << nl;

    dmdt_.writeEntry("dmdt", os);
    dDep_.writeEntry("dDep", os);
    qq_.writeEntry("qQuenching", os);
    alphatConv_.writeEntry("alphatConv", os);
    writeEntry("value", os);
}


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

makePatchTypeField
(
    fvPatchScalarField,
    alphatWallBoilingWallFunctionFvPatchScalarField
);

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace compressible
} // End namespace Foam

// ************************************************************************* //