parcel.C 19.3 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
Mark Olesen's avatar
Mark Olesen committed
5
    \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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
     \\/     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 2 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, write to the Free Software Foundation,
    Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA

\*---------------------------------------------------------------------------*/

#include "parcel.H"

#include "spray.H"
#include "dragModel.H"
#include "evaporationModel.H"
#include "heatTransferModel.H"
#include "wallModel.H"
#include "wallPolyPatch.H"
#include "wedgePolyPatch.H"
#include "processorPolyPatch.H"
#include "combustionMixture.H"

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

41
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43
44
45
namespace Foam
{
defineParticleTypeNameAndDebug(parcel, 0);
defineTemplateTypeNameAndDebug(Cloud<parcel>, 0);
46
}
47
48
49

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

50
Foam::parcel::parcel
51
52
53
(
    const Cloud<parcel>& cloud,
    const vector& position,
54
    const label cellI,
55
56
57
58
59
60
61
62
63
64
65
66
67
68
    const vector& n,
    const scalar d,
    const scalar T,
    const scalar m,
    const scalar y,
    const scalar yDot,
    const scalar ct,
    const scalar ms,
    const scalar tTurb,
    const scalar liquidCore,
    const scalar injector,
    const vector& U,
    const vector& Uturb,
    const scalarField& X,
69
    const List<word>& liquidNames
70
71
)
:
72
73
    Particle<parcel>(cloud, position, cellI),
    liquidComponents_
74
    (
75
        liquidNames
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
    ),
    d_(d),
    T_(T),
    m_(m),
    y_(y),
    yDot_(yDot),
    ct_(ct),
    ms_(ms),
    tTurb_(tTurb),
    liquidCore_(liquidCore),
    injector_(injector),
    U_(U),
    Uturb_(Uturb),
    n_(n),
    X_(X),
    tMom_(GREAT)
{}


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

97
bool Foam::parcel::move(spray& sDB)
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
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
291
292
293
294
295
296
297
298
299
{
    const polyMesh& mesh = cloud().pMesh();
    const polyBoundaryMesh& pbMesh = mesh.boundaryMesh();

    const liquidMixture& fuels = sDB.fuels();

    scalar deltaT = sDB.runTime().deltaT().value();
    label Nf = fuels.components().size();
    label Ns = sDB.composition().Y().size();

    // Calculate the interpolated gas properties at the position of the parcel
    vector Up = sDB.UInterpolator().interpolate(position(), cell())
        + Uturb();
    scalar rhog = sDB.rhoInterpolator().interpolate(position(), cell());
    scalar pg = sDB.pInterpolator().interpolate(position(), cell());
    scalar Tg = sDB.TInterpolator().interpolate(position(), cell());

    scalarField Yfg(Nf, 0.0);

    scalar cpMixture = 0.0;
    for(label i=0; i<Ns; i++)
    {
        const volScalarField& Yi = sDB.composition().Y()[i];
        if (sDB.isLiquidFuel()[i])
        {
            label j = sDB.gasToLiquidIndex()[i];
            scalar Yicelli = Yi[cell()];
            Yfg[j] = Yicelli;
        }
        cpMixture += Yi[cell()]*sDB.gasProperties()[i].Cp(Tg);
    }

    // correct the gaseous temperature for evaporated fuel

    scalar cellV = sDB.mesh().V()[cell()];
    scalar cellMass = rhog*cellV;
    Tg += sDB.shs()[cell()]/(cpMixture*cellMass);
    Tg = max(200.0, Tg);

    scalar tauMomentum = GREAT;
    scalar tauHeatTransfer = GREAT;
    scalarField tauEvaporation(Nf, GREAT);
    scalarField tauBoiling(Nf, GREAT);

    bool keepParcel = true;

    setRelaxationTimes
    (
        cell(),
        tauMomentum,
        tauEvaporation,
        tauHeatTransfer,
        tauBoiling,
        sDB,
        rhog,
        Up,
        Tg,
        pg,
        Yfg,
        m()*fuels.Y(X()),
        deltaT
    );


    // set the end-time for the track
    scalar tEnd = (1.0 - stepFraction())*deltaT;

    // set the maximum time step for this parcel
    scalar dtMax = min
    (
        tEnd,
        min
        (
            tauMomentum,
            min
            (
                1.0e+10*mag(min(tauEvaporation)), // evaporation is not an issue
                min
                (
                    mag(tauHeatTransfer),
                    mag(min(tauBoiling))
                )
            )
        )
    )/sDB.subCycles();

    // prevent the number of subcycles from being too many
    // (10 000 seems high enough)
    dtMax = max(dtMax, 1.0e-4*tEnd);

    bool switchProcessor = false;
    vector planeNormal = vector::zero;
    if (sDB.twoD())
    {
        planeNormal = n() ^ sDB.axisOfSymmetry();
        planeNormal /= mag(planeNormal);
    }

    // move the parcel until there is no 'timeLeft'
    while (keepParcel && tEnd > SMALL && !switchProcessor)
    {
        // set the lagrangian time-step
        scalar dt = min(dtMax, tEnd);

        // remember which cell the parcel is in
        // since this will change if a face is hit
        label celli = cell();
        scalar p = sDB.p()[celli];

        // track parcel to face, or end of trajectory
        if (keepParcel)
        {
            // Track and adjust the time step if the trajectory is not completed
            dt *= trackToFace(position() + dt*U_, sDB);

            // Decrement the end-time acording to how much time the track took
            tEnd -= dt;

            // Set the current time-step fraction.
            stepFraction() = 1.0 - tEnd/deltaT;

            if (onBoundary()) // hit face
            {
#               include "boundaryTreatment.H"
            }
        }

        if (keepParcel && sDB.twoD())
        {
            scalar z = position() & sDB.axisOfSymmetry();
            vector r = position() - z*sDB.axisOfSymmetry();
            if (mag(r) > SMALL)
            {
                correctNormal(sDB.axisOfSymmetry());
            }
        }

        // **** calculate the lagrangian source terms ****
        // First we get the 'old' properties.
        // and then 'update' them to get the 'new'
        // properties.
        // The difference is then added to the source terms.

        scalar oRho = fuels.rho(p, T(), X());
        scalarField oMass(Nf, 0.0);
        scalar oHg = 0.0;
        scalar oTotMass = m();
        scalarField oYf(fuels.Y(X()));

        forAll(oMass, i)
        {
            oMass[i] = m()*oYf[i];
            label j = sDB.liquidToGasIndex()[i];
            oHg += oYf[i]*sDB.gasProperties()[j].H(T());
        }

        vector oMom = m()*U();
        scalar oHv = fuels.hl(p, T(), X());
        scalar oH = oHg - oHv;
        scalar oPE = (p - fuels.pv(p, T(), X()))/oRho;

        // update the parcel properties (U, T, D)
        updateParcelProperties
        (
            dt,
            sDB,
            celli,
            face()
        );

        scalar nRho = fuels.rho(p, T(), X());
        scalar nHg = 0.0;
        scalarField nMass(Nf, 0.0);
        scalarField nYf(fuels.Y(X()));

        forAll(nMass, i)
        {
            nMass[i] = m()*nYf[i];
            label j = sDB.liquidToGasIndex()[i];
            nHg += nYf[i]*sDB.gasProperties()[j].H(T());
        }

        vector nMom = m()*U();
        scalar nHv = fuels.hl(p, T(), X());
        scalar nH = nHg - nHv;
        scalar nPE = (p - fuels.pv(p, T(), X()))/nRho;

        // Update the Spray Source Terms
        forAll(nMass, i)
        {
            sDB.srhos()[i][celli] += oMass[i] - nMass[i];
        }
        sDB.sms()[celli]   += oMom - nMom;

        sDB.shs()[celli]   +=
            oTotMass*(oH + oPE)
          - m()*(nH + nPE);

        // Remove evaporated mass from stripped mass
        ms() -= ms()*(oTotMass-m())/oTotMass;

        // remove parcel if it is 'small'
henry's avatar
henry committed
300
        if (m() < 1.0e-12)
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
        {
            keepParcel = false;

            // ... and add the removed 'stuff' to the gas
            forAll(nMass, i)
            {
                sDB.srhos()[i][celli] += nMass[i];
            }

            sDB.sms()[celli] += nMom;
            sDB.shs()[celli] += m()*(nH + nPE);
        }

        if (onBoundary() && keepParcel)
        {
            if (face() > -1)
            {
                if (isType<processorPolyPatch>(pbMesh[patch(face())]))
                {
                    switchProcessor = true;
                }
            }
        }
    }

    return keepParcel;
}


330
void Foam::parcel::updateParcelProperties
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
(
    const scalar dt,
    spray& sDB,
    const label celli,
    const label facei
)
{
    const liquidMixture& fuels = sDB.fuels();

    label Nf = sDB.fuels().components().size();
    label Ns = sDB.composition().Y().size();

    // calculate mean molecular weight
    scalar W = 0.0;
    for(label i=0; i<Ns; i++)
    {
        W += sDB.composition().Y()[i][celli]/sDB.gasProperties()[i].W();

    }
    W = 1.0/W;

    // Calculate the interpolated gas properties at the position of the parcel
    vector Up = sDB.UInterpolator().interpolate(position(), celli, facei)
        + Uturb();
    scalar rhog = sDB.rhoInterpolator().interpolate(position(), celli, facei);
    scalar pg = sDB.pInterpolator().interpolate(position(), celli, facei);
    scalar Tg0 = sDB.TInterpolator().interpolate(position(), celli, facei);

    // correct the gaseous temperature for evaporated fuel
    scalar cpMix = 0.0;
    for(label i=0; i<Ns; i++)
    {
        cpMix += sDB.composition().Y()[i][celli]
                *sDB.gasProperties()[i].Cp(T());
    }
    scalar cellV            = sDB.mesh().V()[celli];
    scalar rho              = sDB.rho()[celli];
    scalar cellMass         = rho*cellV;
    scalar dh               = sDB.shs()[celli];
    scalarField addedMass(Nf, 0.0);

    forAll(addedMass, i)
    {
        addedMass[i] += sDB.srhos()[i][celli]*cellV;
    }

    scalar Tg  = Tg0 + dh/(cpMix*cellMass);
    Tg = max(200.0, Tg);

    scalarField Yfg(Nf, 0.0);
    forAll(Yfg, i)
    {
        label j = sDB.liquidToGasIndex()[i];
        const volScalarField& Yj = sDB.composition().Y()[j];
        scalar Yfg0 = Yj[celli];
        Yfg[i] = (Yfg0*cellMass + addedMass[i])/(addedMass[i] + cellMass);
    }

    scalar tauMomentum     = GREAT;
    scalar tauHeatTransfer = GREAT;
    scalarField tauEvaporation(Nf, GREAT);
    scalarField tauBoiling(Nf, GREAT);

    setRelaxationTimes
    (
        celli,
        tauMomentum,
        tauEvaporation,
        tauHeatTransfer,
        tauBoiling,
        sDB,
        rhog,
        Up,
        Tg,
        pg,
        Yfg,
        m()*fuels.Y(X()),
        dt
    );

    scalar timeRatio = dt/tauMomentum;

    vector Ucorr = Up;
    vector gcorr = sDB.g();

    if (sDB.twoD())
    {
        // remove the tangential velocity component
        scalar v1 = Up & sDB.axisOfSymmetry();
        scalar v2 = Up & n();
        Ucorr     = v1*sDB.axisOfSymmetry() + v2*n();

        // Remove the tangential gravity component
        scalar g1 = gcorr & sDB.axisOfSymmetry();
        scalar g2 = gcorr & n();
        gcorr     = g1*sDB.axisOfSymmetry() + g2*n();
    }

    U() = (U() + timeRatio*Ucorr + gcorr*dt)/(1.0 + timeRatio);

    if (sDB.twoD())
    {
        vector normal = n() ^ sDB.axisOfSymmetry();
        normal /= mag(normal);
        scalar dU = U() & normal;
        U() -= dU*normal;
    }

    scalar TDroplet = T();
    scalar oldDensity = fuels.rho(pg, T(), X());
    scalar oldMass = m();
    scalarField Yf0(fuels.Y(X()));
    scalarField mi0(m()*Yf0);
    scalarField mi(mi0);

    scalar oldhg = 0.0;
    for (label i=0; i<Nf; i++)
    {
        label j = sDB.liquidToGasIndex()[i];
        oldhg += Yf0[i]*sDB.gasProperties()[j].H(T());
    }

    scalar oldhv = fuels.hl(pg, T(), X());
    scalar Np = N(oldDensity);

    scalar newDensity = oldDensity;
    scalar newMass    = oldMass;
    scalar newhg      = oldhg;
    scalar newhv      = oldhv;

    scalar Tnew = T();

    // NN.
    // first calculate the new temperature and droplet mass,
    // then calculate the energy source and correct the
    // gaseous temperature, Tg, and mass fraction, Yfg,
    // to calculate the new properties for the parcel
    // This procedure seems to be more stable
    label n = 0;
    while ((n < sDB.evaporation().nEvapIter()) && (m() > VSMALL))
    {
        n++;
        // new characteristic times does not need to be calculated the first time
        if (n > 1)
        {
            newDensity = fuels.rho(pg, Tnew, X());
            newMass = m();
            newhg = 0.0;
            scalarField Ynew(fuels.Y(X()));
            for(label i=0; i<Nf; i++)
            {
                label j = sDB.liquidToGasIndex()[i];
                newhg += Ynew[i]*sDB.gasProperties()[j].H(Tnew);
            }

            newhv = fuels.hl(pg, Tnew, X());

            scalar dm = oldMass - newMass;
            scalar dhNew = oldMass*(oldhg-oldhv) - newMass*(newhg-newhv);

            // Prediction of new gaseous temperature and fuel mass fraction
            Tg  = Tg0 + (dh+dhNew)/(cpMix*cellMass);
            Tg = max(200.0, Tg);

            forAll(Yfg, i)
            {
                label j = sDB.liquidToGasIndex()[i];
                const volScalarField& Yj = sDB.composition().Y()[j];
                scalar Yfg0 = Yj[celli];
                Yfg[i] = (Yfg0*cellMass + addedMass[i] + dm)
                        /(addedMass[i] + cellMass + dm);
            }

            setRelaxationTimes
            (
                celli,
                tauMomentum,
                tauEvaporation,
                tauHeatTransfer,
                tauBoiling,
                sDB,
                rhog,
                Up,
                Tg,
                pg,
                Yfg,
                m()*fuels.Y(X()),
                dt
            );

        }

        scalar Taverage = TDroplet + (Tg - TDroplet)/3.0;
        // for a liquid Cl \approx Cp
        scalar liquidcL = sDB.fuels().cp(pg, TDroplet, X());

        cpMix = 0.0;
        for(label i=0; i<Ns; i++)
        {
            if (sDB.isLiquidFuel()[i])
            {
                label j = sDB.gasToLiquidIndex()[i];
                cpMix += Yfg[j]*sDB.gasProperties()[i].Cp(Taverage);
            }
            else
            {
                scalar Y = sDB.composition().Y()[i][celli];
                cpMix += Y*sDB.gasProperties()[i].Cp(Taverage);
            }
        }

        scalar evaporationSource = 0.0;
        scalar z = 0.0;
        scalar tauEvap = 0.0;

        for(label i=0; i<Nf; i++)
        {
            tauEvap += X()[i]*fuels.properties()[i].W()/tauEvaporation[i];
        }
        tauEvap = fuels.W(X())/tauEvap;


        if (sDB.evaporation().evaporation())
        {
            scalar hv = fuels.hl(pg, TDroplet, X());
            evaporationSource =
                hv/liquidcL/tauEvap;

            z = cpMix*tauHeatTransfer/liquidcL/tauEvap;
        }

        if (sDB.heatTransfer().heatTransfer())
        {
            scalar fCorrect =
                sDB.heatTransfer().fCorrection(z)/tauHeatTransfer;

            Tnew =
                (TDroplet + dt*(fCorrect * Tg - evaporationSource))
                /(1.0 + dt*fCorrect);

            // Prevent droplet temperature to go above critial value
            Tnew = min(Tnew, fuels.Tc(X()));

            // Prevent droplet temperature to go too low
            // Mainly a numerical stability issue
            Tnew = max(200.0, Tnew);
henry's avatar
henry committed
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
            scalar Td = Tnew;

            scalar pAtSurface = fuels.pv(pg, Td, X());
            scalar pCompare = 0.999*pg;
            scalar boiling = pAtSurface >= pCompare;
            if (boiling)
            {
                // can not go above boiling temperature
                scalar Terr = 1.0e-3;
                label n=0;
                scalar dT = 1.0;
                scalar pOld = pAtSurface;
                while (dT > Terr)
                {
                    n++;
                    pAtSurface = fuels.pv(pg, Td, X());
                    if ((pAtSurface < pCompare) && (pOld < pCompare))
                    {
                        Td += dT;
                    }
                    else
                    {
                        if ((pAtSurface > pCompare) && (pOld > pCompare))
                        {
                            Td -= dT;
                        }
                        else
                        {
                            dT *= 0.5;
                            if ((pAtSurface > pCompare) && (pOld < pCompare))
                            {
                                Td -= dT;
                            }
                            else
                            {
                                Td += dT;
                            }
                        }
                    }
                    pOld = pAtSurface;
                    if (debug)
                    {
                        if (n>100)
                        {
                            Info << "n = " << n << ", T = " << Td << ", pv = " << pAtSurface << endl;
                        }
                    }
                }
                Tnew = Td;
            }
627
628
629
630
631
632
633
634
635
636
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
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
        }

        // Evaporate droplet!
        // if the droplet is NOT boiling use implicit scheme.
        if (sDB.evaporation().evaporation())
        {
            for(label i=0; i<Nf; i++)
            {
                // immediately evaporate mass that has reached critical
                // condition
                if (mag(Tnew - fuels.Tc(X())) < SMALL)
                {
                    mi[i] = 0.0;
                }
                else
                {
                    scalar Td = min(Tnew, 0.999*fuels.properties()[i].Tc());

                    scalar pAtSurface = fuels.properties()[i].pv(pg, Td);
                    scalar boiling = pAtSurface >= 0.999*pg;

                    if (!boiling)
                    {
                        scalar fr = dt/tauEvaporation[i];
                        mi[i] = mi0[i]/(1.0 + fr);
                    }
                    else
                    {
                        scalar fr = dt/tauBoiling[i];
                        mi[i] = mi0[i]/(1.0 + fr);
                    }
                }
            }

            scalar mTot = sum(mi);
            if (mTot > VSMALL)
            {
                scalarField Ynew(mi/mTot);
                scalarField Xnew(sDB.fuels().X(Ynew));
                forAll(Xnew, i)
                {
                    X()[i] = Xnew[i];
                }
                m() = mTot;
            }
            else
            {
                m() = 0.0;
            }
        }
        T() = Tnew;
        scalar rhod = fuels.rho(pg, T(), X());
        d() = cbrt(6.0*m()/(Np*rhod*M_PI));
    }

    T() = Tnew;
    scalar rhod = fuels.rho(pg, T(), X());
    m() = sum(mi);
    d() = cbrt(6.0*m()/(Np*rhod*M_PI));
}


689
void Foam::parcel::transformProperties(const tensor& T)
690
691
692
693
694
{
    U_ = transform(T, U_);
}


695
void Foam::parcel::transformProperties(const vector&)
696
697
698
699
{}


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