Skip to content
Snippets Groups Projects
epsilonLowReWallFunctionFvPatchScalarField.C 5.1 KiB
Newer Older
  • Learn to ignore specific revisions
  • Henry's avatar
    Henry committed
    /*---------------------------------------------------------------------------*\
      =========                 |
      \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
       \\    /   O peration     |
    
    Henry's avatar
    Henry committed
        \\  /    A nd           | Copyright (C) 2012-2015 OpenFOAM Foundation
    
    Henry's avatar
    Henry committed
         \\/     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 "epsilonLowReWallFunctionFvPatchScalarField.H"
    #include "turbulenceModel.H"
    #include "fvPatchFieldMapper.H"
    #include "volFields.H"
    #include "addToRunTimeSelectionTable.H"
    
    // * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //
    
    
    Foam::scalar Foam::epsilonLowReWallFunctionFvPatchScalarField::yPlusLam
    
    Henry's avatar
    Henry committed
    (
        const scalar kappa,
        const scalar E
    )
    {
        scalar ypl = 11.0;
    
        for (int i=0; i<10; i++)
        {
            ypl = log(max(E*ypl, 1))/kappa;
        }
    
        return ypl;
    }
    
    
    
    void Foam::epsilonLowReWallFunctionFvPatchScalarField::calculate
    
    Henry's avatar
    Henry committed
    (
    
    Henry's avatar
    Henry committed
        const List<scalar>& cornerWeights,
        const fvPatch& patch,
    
    Henry's avatar
    Henry committed
    )
    {
        const label patchi = patch.index();
    
    
    Henry's avatar
    Henry committed
    
        const scalar Cmu25 = pow025(Cmu_);
        const scalar Cmu75 = pow(Cmu_, 0.75);
    
    
    Henry's avatar
    Henry committed
        const volScalarField& k = tk();
    
    
        const tmp<scalarField> tnuw = turbModel.nu(patchi);
    
    Henry's avatar
    Henry committed
        const scalarField& nuw = tnuw();
    
    
        const tmp<scalarField> tnutw = turbModel.nut(patchi);
    
    Henry's avatar
    Henry committed
        const scalarField& nutw = tnutw();
    
    
        const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
    
    Henry's avatar
    Henry committed
    
        const scalarField magGradUw(mag(Uw.snGrad()));
    
    
        typedef DimensionedField<scalar, volMesh> FieldType;
    
        const DimensionedField<scalar, volMesh>& G =
            db().lookupObject<DimensionedField<scalar, volMesh> >
            (
                turbModel.GName()
            );
    
    
    Henry's avatar
    Henry committed
        // Set epsilon and G
    
    Henry's avatar
    Henry committed
        forAll(nutw, facei)
    
    Henry's avatar
    Henry committed
        {
    
    Henry's avatar
    Henry committed
            label celli = patch.faceCells()[facei];
    
    Henry's avatar
    Henry committed
    
    
    Henry's avatar
    Henry committed
            scalar yPlus = Cmu25*sqrt(k[celli])*y[facei]/nuw[facei];
    
    Henry's avatar
    Henry committed
    
    
    Henry's avatar
    Henry committed
            scalar w = cornerWeights[facei];
    
    Henry's avatar
    Henry committed
    
            if (yPlus > yPlusLam_)
            {
    
                epsilon0[celli] += w*Cmu75*pow(k[celli], 1.5)/(kappa_*y[facei]);
    
                G0[celli] +=
                    w
                   *(nutw[facei] + nuw[facei])
                   *magGradUw[facei]
                   *Cmu25*sqrt(k[celli])
                   /(kappa_*y[facei]);
    
    Henry's avatar
    Henry committed
            }
            else
            {
    
                epsilon0[celli] += w*2.0*k[celli]*nuw[facei]/sqr(y[facei]);
                G0[celli] += G[celli];
    
    Henry's avatar
    Henry committed
            }
        }
    }
    
    
    // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
    
    
    Foam::epsilonLowReWallFunctionFvPatchScalarField::
    
    Henry's avatar
    Henry committed
    epsilonLowReWallFunctionFvPatchScalarField
    (
        const fvPatch& p,
        const DimensionedField<scalar, volMesh>& iF
    )
    :
        epsilonWallFunctionFvPatchScalarField(p, iF),
        yPlusLam_(yPlusLam(kappa_, E_))
    {}
    
    
    
    Foam::epsilonLowReWallFunctionFvPatchScalarField::
    
    Henry's avatar
    Henry committed
    epsilonLowReWallFunctionFvPatchScalarField
    (
        const epsilonLowReWallFunctionFvPatchScalarField& ptf,
        const fvPatch& p,
        const DimensionedField<scalar, volMesh>& iF,
        const fvPatchFieldMapper& mapper
    )
    :
        epsilonWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
        yPlusLam_(ptf.yPlusLam_)
    {}
    
    
    
    Foam::epsilonLowReWallFunctionFvPatchScalarField::
    
    Henry's avatar
    Henry committed
    epsilonLowReWallFunctionFvPatchScalarField
    (
        const fvPatch& p,
        const DimensionedField<scalar, volMesh>& iF,
        const dictionary& dict
    )
    :
        epsilonWallFunctionFvPatchScalarField(p, iF, dict),
        yPlusLam_(yPlusLam(kappa_, E_))
    {}
    
    
    
    Foam::epsilonLowReWallFunctionFvPatchScalarField::
    
    Henry's avatar
    Henry committed
    epsilonLowReWallFunctionFvPatchScalarField
    (
        const epsilonLowReWallFunctionFvPatchScalarField& ewfpsf
    )
    :
        epsilonWallFunctionFvPatchScalarField(ewfpsf),
        yPlusLam_(ewfpsf.yPlusLam_)
    {}
    
    
    
    Foam::epsilonLowReWallFunctionFvPatchScalarField::
    
    Henry's avatar
    Henry committed
    epsilonLowReWallFunctionFvPatchScalarField
    (
        const epsilonLowReWallFunctionFvPatchScalarField& ewfpsf,
        const DimensionedField<scalar, volMesh>& iF
    )
    :
        epsilonWallFunctionFvPatchScalarField(ewfpsf, iF),
        yPlusLam_(ewfpsf.yPlusLam_)
    {}
    
    
    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
    
    
    namespace Foam
    {
        makePatchTypeField
        (
            fvPatchScalarField,
            epsilonLowReWallFunctionFvPatchScalarField
        );
    }
    
    Henry's avatar
    Henry committed
    
    
    // ************************************************************************* //