alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C 10 KB
Newer Older
1 2 3 4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2015-2016 OpenFOAM Foundation
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
     \\/     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 "alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H"
#include "fvPatchFieldMapper.H"
#include "addToRunTimeSelectionTable.H"

#include "twoPhaseSystem.H"
#include "ThermalPhaseChangePhaseSystem.H"
#include "MomentumTransferPhaseSystem.H"
#include "compressibleTurbulenceModel.H"
#include "ThermalDiffusivity.H"
#include "PhaseCompressibleTurbulenceModel.H"
#include "wallFvPatch.H"

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

namespace Foam
{
namespace compressible
{

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

scalar alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxExp_
    = 50.0;
scalar alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::tolerance_
    = 0.01;
label alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxIters_
    = 10;

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

void alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::checkType()
{
    if (!isA<wallFvPatch>(patch()))
    {
        FatalErrorInFunction
            << "Patch type for patch " << patch().name() << " must be wall\n"
            << "Current patch type is " << patch().type() << nl
            << exit(FatalError);
    }
}


tmp<scalarField>
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Psmooth
(
    const scalarField& Prat
) const
{
    return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
}


tmp<scalarField>
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
(
    const scalarField& P,
    const scalarField& Prat
) const
{
    tmp<scalarField> typsf(new scalarField(this->size()));
86
    scalarField& ypsf = typsf.ref();
87

88
    forAll(ypsf, facei)
89 90 91 92 93
    {
        scalar ypt = 11.0;

        for (int i=0; i<maxIters_; i++)
        {
94 95
            scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
            scalar df = 1 - 1.0/(ypt*kappa_*Prat[facei]);
96 97 98 99
            scalar yptNew = ypt - f/df;

            if (yptNew < VSMALL)
            {
100
                ypsf[facei] = 0;
101 102 103
            }
            else if (mag(yptNew - ypt) < tolerance_)
            {
104
                ypsf[facei] = yptNew;
105 106 107 108 109 110 111
            }
            else
            {
                ypt = yptNew;
            }
        }

112
        ypsf[facei] = ypt;
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
    }

    return typsf;
}

tmp<scalarField>
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::calcAlphat
(
    const scalarField& prevAlphat
) const
{
    // Lookup the fluid model
    const ThermalPhaseChangePhaseSystem
    <
        MomentumTransferPhaseSystem<twoPhaseSystem>
    >& fluid =
        refCast
        <
            const ThermalPhaseChangePhaseSystem
            <
                MomentumTransferPhaseSystem<twoPhaseSystem>
            >
        >
        (
            db().lookupObject<phaseSystem>("phaseProperties")
        );

    const phaseModel& liquid
    (
142
        fluid.phase1().name() == internalField().group()
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
      ? fluid.phase1()
      : fluid.phase2()
    );

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

    // Retrieve turbulence properties from model
    const phaseCompressibleTurbulenceModel& turbModel = liquid.turbulence();

    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().alpha(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& hew =
        liquid.thermo().he().boundaryField()[patchi];

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

    scalarField Tp(Tw.patchInternalField());

    // Heat flux [W/m2] - lagging alphatw
    const scalarField qDot
    (
        (prevAlphat + alphaw)*hew.snGrad()
    );

    scalarField uTau(Cmu25*sqrt(kw));

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

    scalarField Pr(muw/alphaw);

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

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

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

    tmp<scalarField> talphatConv(new scalarField(this->size()));
200
    scalarField& alphatConv = talphatConv.ref();
201 202

    // Populate boundary values
203
    forAll(alphatConv, facei)
204 205 206
    {
        // Evaluate new effective thermal diffusivity
        scalar alphaEff = 0.0;
207
        if (yPlus[facei] < yPlusTherm[facei])
208
        {
209 210 211
            scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
            scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
            scalar C = Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
212 213 214 215
            alphaEff = A/(B + C + VSMALL);
        }
        else
        {
216
            scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
217
            scalar B =
218
                qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
219
            scalar magUc =
220
                uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
221
            scalar C =
222 223
                0.5*rhow[facei]*uTau[facei]
               *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
224 225 226 227
            alphaEff = A/(B + C + VSMALL);
        }

        // Update convective heat transfer turbulent thermal diffusivity
228
        alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
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 300 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 330 331 332 333 334 335 336
    }

    return talphatConv;
}


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

alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
(
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF
)
:
    alphatPhaseChangeWallFunctionFvPatchScalarField(p, iF),
    Prt_(0.85),
    Cmu_(0.09),
    kappa_(0.41),
    E_(9.8)
{
    checkType();
}


alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
(
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF,
    const dictionary& dict
)
:
    alphatPhaseChangeWallFunctionFvPatchScalarField(p, iF, dict),
    Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
    Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
    kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
    E_(dict.lookupOrDefault<scalar>("E", 9.8))
{}


alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
(
    const alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField& ptf,
    const fvPatch& p,
    const DimensionedField<scalar, volMesh>& iF,
    const fvPatchFieldMapper& mapper
)
:
    alphatPhaseChangeWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
    Prt_(ptf.Prt_),
    Cmu_(ptf.Cmu_),
    kappa_(ptf.kappa_),
    E_(ptf.E_)
{}


alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
(
    const alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField& awfpsf
)
:
    alphatPhaseChangeWallFunctionFvPatchScalarField(awfpsf),
    Prt_(awfpsf.Prt_),
    Cmu_(awfpsf.Cmu_),
    kappa_(awfpsf.kappa_),
    E_(awfpsf.E_)
{}


alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
(
    const alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField& awfpsf,
    const DimensionedField<scalar, volMesh>& iF
)
:
    alphatPhaseChangeWallFunctionFvPatchScalarField(awfpsf, iF),
    Prt_(awfpsf.Prt_),
    Cmu_(awfpsf.Cmu_),
    kappa_(awfpsf.kappa_),
    E_(awfpsf.E_)
{}


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

void alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::updateCoeffs()
{
    if (updated())
    {
        return;
    }

    operator==(calcAlphat(*this));

    fixedValueFvPatchScalarField::updateCoeffs();
}


void alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::write
(
    Ostream& os
) const
{
    fvPatchField<scalar>::write(os);
337 338 339 340
    os.writeEntry("Prt", Prt_);
    os.writeEntry("Cmu", Cmu_);
    os.writeEntry("kappa", kappa_);
    os.writeEntry("E", E_);
341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360
    dmdt_.writeEntry("dmdt", os);
    writeEntry("value", os);
}


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

makePatchTypeField
(
    fvPatchScalarField,
    alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
);


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

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

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