ReynoldsAnalogy.C 7.35 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
     \\/     M anipulation  |
OpenFOAM bot's avatar
OpenFOAM bot committed
7
-------------------------------------------------------------------------------
8
    Copyright (C) 2017-2020 OpenCFD Ltd.
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
-------------------------------------------------------------------------------
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 "ReynoldsAnalogy.H"
#include "fluidThermo.H"
#include "turbulentTransportModel.H"
#include "turbulentFluidThermoModel.H"
#include "addToRunTimeSelectionTable.H"

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

namespace Foam
{
namespace heatTransferCoeffModels
{
    defineTypeNameAndDebug(ReynoldsAnalogy, 0);
    addToRunTimeSelectionTable
    (
        heatTransferCoeffModel,
        ReynoldsAnalogy,
        dictionary
    );
}
}


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

Foam::tmp<Foam::Field<Foam::scalar>>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::rho(const label patchi) const
{
    if (rhoName_ == "rhoInf")
    {
        const label n = mesh_.boundary()[patchi].size();
59
        return tmp<Field<scalar>>::New(n, rhoRef_);
60 61 62 63 64 65 66 67
    }
    else if (mesh_.foundObject<volScalarField>(rhoName_, false))
    {
        const volScalarField& rho =
            mesh_.lookupObject<volScalarField>(rhoName_);
        return rho.boundaryField()[patchi];
    }

68 69 70 71
    FatalErrorInFunction
        << "Unable to set rho for patch " << patchi
        << exit(FatalError);

72
    return nullptr;
73 74 75 76 77 78 79 80 81
}


Foam::tmp<Foam::Field<Foam::scalar>>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::Cp(const label patchi) const
{
    if (CpName_ == "CpInf")
    {
        const label n = mesh_.boundary()[patchi].size();
82
        return tmp<Field<scalar>>::New(n, CpRef_);
83
    }
84
    else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
85 86
    {
        const fluidThermo& thermo =
87
            mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);
88 89 90 91 92 93 94

        const scalarField& pp = thermo.p().boundaryField()[patchi];
        const scalarField& Tp = thermo.T().boundaryField()[patchi];

        return thermo.Cp(pp, Tp, patchi);
    }

95 96 97 98
    FatalErrorInFunction
        << "Unable to set Cp for patch " << patchi
        << exit(FatalError);

99
    return nullptr;
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
}


Foam::tmp<Foam::volSymmTensorField>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::devReff() const
{
    typedef compressible::turbulenceModel cmpTurbModel;
    typedef incompressible::turbulenceModel icoTurbModel;

    if (mesh_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
    {
        const cmpTurbModel& turb =
            mesh_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);

        return turb.devRhoReff()/turb.rho();
    }
    else if (mesh_.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
    {
        const incompressible::turbulenceModel& turb =
            mesh_.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);

        return turb.devReff();
    }
    else if (mesh_.foundObject<fluidThermo>(fluidThermo::dictName))
    {
        const fluidThermo& thermo =
            mesh_.lookupObject<fluidThermo>(fluidThermo::dictName);

        const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);

        return -thermo.nu()*dev(twoSymm(fvc::grad(U)));
    }
    else if (mesh_.foundObject<transportModel>("transportProperties"))
    {
        const transportModel& laminarT =
            mesh_.lookupObject<transportModel>("transportProperties");

        const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);

        return -laminarT.nu()*dev(twoSymm(fvc::grad(U)));
    }
    else if (mesh_.foundObject<dictionary>("transportProperties"))
    {
        const dictionary& transportProperties =
            mesh_.lookupObject<dictionary>("transportProperties");

146
        dimensionedScalar nu("nu", dimViscosity, transportProperties);
147 148 149 150 151 152

        const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);

        return -nu*dev(twoSymm(fvc::grad(U)));
    }

153 154 155 156 157
    FatalErrorInFunction
        << "No valid model for viscous stress calculation"
        << exit(FatalError);

    return nullptr;
158 159 160 161 162 163 164 165 166
}


Foam::tmp<Foam::FieldField<Foam::Field, Foam::scalar>>
Foam::heatTransferCoeffModels::ReynoldsAnalogy::Cf() const
{
    const volVectorField& U = mesh_.lookupObject<volVectorField>(UName_);
    const volVectorField::Boundary& Ubf = U.boundaryField();

167 168
    auto tCf = tmp<FieldField<Field, scalar>>::New(Ubf.size());
    auto& Cf = tCf.ref();
169 170 171

    forAll(Cf, patchi)
    {
172
        Cf.set(patchi, new Field<scalar>(Ubf[patchi].size(), Zero));
173 174 175 176 177
    }

    const volSymmTensorField R(devReff());
    const volSymmTensorField::Boundary& Rbf = R.boundaryField();

178
    for (const label patchi : patchSet_)
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
    {
        const fvPatchVectorField& Up = Ubf[patchi];

        const symmTensorField& Rp = Rbf[patchi];

        const vectorField nHat(Up.patch().nf());

        const scalarField tauByRhop(mag(nHat & Rp));

        Cf[patchi] = 2*tauByRhop/magSqr(URef_);
    }

    return tCf;
}


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

Foam::heatTransferCoeffModels::ReynoldsAnalogy::ReynoldsAnalogy
(
    const dictionary& dict,
    const fvMesh& mesh,
    const word& TName
)
:
    heatTransferCoeffModel(dict, mesh, TName),
    UName_("U"),
206
    URef_(Zero),
207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224
    rhoName_("rho"),
    rhoRef_(0.0),
    CpName_("Cp"),
    CpRef_(0.0)
{
    read(dict);
}


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

bool Foam::heatTransferCoeffModels::ReynoldsAnalogy::read
(
    const dictionary& dict
)
{
    if (heatTransferCoeffModel::read(dict))
    {
225
        dict.readEntry("UInf", URef_);
226 227 228 229

        dict.readIfPresent("Cp", CpName_);
        if (CpName_ == "CpInf")
        {
230
            dict.readEntry("CpInf", CpRef_);
231 232 233 234 235
        }

        dict.readIfPresent("rho", rhoName_);
        if (rhoName_ == "rhoInf")
        {
236
            dict.readEntry("rhoInf", rhoRef_);
237 238 239 240 241 242 243 244 245
        }

        return true;
    }

    return false;
}


246 247 248 249 250
void Foam::heatTransferCoeffModels::ReynoldsAnalogy::htc
(
    volScalarField& htc,
    const FieldField<Field, scalar>& q
)
251 252 253 254 255
{
    const FieldField<Field, scalar> CfBf(Cf());
    const scalar magU = mag(URef_);

    volScalarField::Boundary& htcBf = htc.boundaryFieldRef();
256 257

    for (const label patchi : patchSet_)
258 259 260 261 262 263 264 265 266 267
    {
        const scalarField rhop(rho(patchi));
        const scalarField Cpp(Cp(patchi));

        htcBf[patchi] = 0.5*rhop*Cpp*magU*CfBf[patchi];
    }
}


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