Added mutURoughWallFunction for compatibility with incompressible wall functions
This BC is the same as mutSpalartAlmarasStandardRoughWallFunction. We will refactor mutSpalartAlmarasStandardRoughWallFunction such that it derives from this BC in order to have both backward compatibility and compatibility with incompressible wall functions.
This commit is contained in:
parent
660326047d
commit
6681c7e1ee
3 changed files with 554 additions and 0 deletions
|
@ -22,6 +22,7 @@ $(mutWallFunctions)/mutLowReWallFunction/mutLowReWallFunctionFvPatchScalarField.
|
|||
$(mutWallFunctions)/mutkWallFunction/mutkWallFunctionFvPatchScalarField.C
|
||||
$(mutWallFunctions)/mutkRoughWallFunction/mutkRoughWallFunctionFvPatchScalarField.C
|
||||
$(mutWallFunctions)/mutUWallFunction/mutUWallFunctionFvPatchScalarField.C
|
||||
$(mutWallFunctions)/mutURoughWallFunction/mutURoughWallFunctionFvPatchScalarField.C
|
||||
|
||||
|
||||
epsilonWallFunctions = $(wallFunctions)/epsilonWallFunctions
|
||||
|
|
|
@ -0,0 +1,316 @@
|
|||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | foam-extend: Open Source CFD
|
||||
\\ / O peration | Version: 4.1
|
||||
\\ / A nd | Web: http://www.foam-extend.org
|
||||
\\/ M anipulation | For copyright notice see file Copyright
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of foam-extend.
|
||||
|
||||
foam-extend 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.
|
||||
|
||||
foam-extend 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 foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "mutURoughWallFunctionFvPatchScalarField.H"
|
||||
#include "RASModel.H"
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace compressible
|
||||
{
|
||||
namespace RASModels
|
||||
{
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
tmp<scalarField> mutURoughWallFunctionFvPatchScalarField::calcYPlus
|
||||
(
|
||||
const scalarField& magUp
|
||||
) const
|
||||
{
|
||||
const label patchI = patch().index();
|
||||
|
||||
const turbulenceModel& turbModel =
|
||||
db().lookupObject<turbulenceModel>("turbulenceModel");
|
||||
|
||||
const scalarField& y = turbModel.y()[patchI];
|
||||
const scalarField& muw = turbModel.mu().boundaryField()[patchI];
|
||||
const fvPatchScalarField& rho = turbModel.rho().boundaryField()[patchI];
|
||||
|
||||
tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
|
||||
scalarField& yPlus = tyPlus();
|
||||
|
||||
if (roughnessHeight_ > 0.0)
|
||||
{
|
||||
// Rough Walls
|
||||
const scalar c_1 = 1/(90 - 2.25) + roughnessConstant_;
|
||||
static const scalar c_2 = 2.25/(90 - 2.25);
|
||||
static const scalar c_3 = 2.0*atan(1.0)/log(90/2.25);
|
||||
static const scalar c_4 = c_3*log(2.25);
|
||||
|
||||
//if (KsPlusBasedOnYPlus_)
|
||||
{
|
||||
// If KsPlus is based on YPlus the extra term added to the law
|
||||
// of the wall will depend on yPlus
|
||||
forAll(yPlus, faceI)
|
||||
{
|
||||
const scalar magUpara = magUp[faceI];
|
||||
const scalar Re = rho[faceI]*magUpara*y[faceI]/muw[faceI];
|
||||
const scalar kappaRe = kappa_*Re;
|
||||
|
||||
scalar yp = yPlusLam_;
|
||||
const scalar ryPlusLam = 1.0/yp;
|
||||
|
||||
label iter = 0;
|
||||
scalar yPlusLast = 0.0;
|
||||
scalar dKsPlusdYPlus = roughnessHeight_/y[faceI];
|
||||
|
||||
// Enforce the roughnessHeight to be less than the distance to
|
||||
// the first cell centre.
|
||||
if (dKsPlusdYPlus > 1)
|
||||
{
|
||||
dKsPlusdYPlus = 1;
|
||||
}
|
||||
|
||||
// Additional tuning parameter (fudge factor) - nominally = 1
|
||||
dKsPlusdYPlus *= roughnessFactor_;
|
||||
|
||||
do
|
||||
{
|
||||
yPlusLast = yp;
|
||||
|
||||
// The non-dimensional roughness height
|
||||
scalar KsPlus = yp*dKsPlusdYPlus;
|
||||
|
||||
// The extra term in the law-of-the-wall
|
||||
scalar G = 0.0;
|
||||
|
||||
scalar yPlusGPrime = 0.0;
|
||||
|
||||
if (KsPlus >= 90)
|
||||
{
|
||||
const scalar t_1 = 1 + roughnessConstant_*KsPlus;
|
||||
G = log(t_1);
|
||||
yPlusGPrime = roughnessConstant_*KsPlus/t_1;
|
||||
}
|
||||
else if (KsPlus > 2.25)
|
||||
{
|
||||
const scalar t_1 = c_1*KsPlus - c_2;
|
||||
const scalar t_2 = c_3*log(KsPlus) - c_4;
|
||||
const scalar sint_2 = sin(t_2);
|
||||
const scalar logt_1 = log(t_1);
|
||||
G = logt_1*sint_2;
|
||||
yPlusGPrime =
|
||||
(c_1*sint_2*KsPlus/t_1) + (c_3*logt_1*cos(t_2));
|
||||
}
|
||||
|
||||
scalar denom = 1.0 + log(E_*yp) - G - yPlusGPrime;
|
||||
if (mag(denom) > VSMALL)
|
||||
{
|
||||
yp = (kappaRe + yp*(1 - yPlusGPrime))/denom;
|
||||
}
|
||||
} while
|
||||
(
|
||||
mag(ryPlusLam*(yp - yPlusLast)) > 0.0001
|
||||
&& ++iter < 10
|
||||
&& yp > VSMALL
|
||||
);
|
||||
|
||||
yPlus[faceI] = max(0.0, yp);
|
||||
}
|
||||
}
|
||||
}
|
||||
else
|
||||
{
|
||||
// Smooth Walls
|
||||
forAll(yPlus, faceI)
|
||||
{
|
||||
const scalar magUpara = magUp[faceI];
|
||||
const scalar Re = rho[faceI]*magUpara*y[faceI]/muw[faceI];
|
||||
const scalar kappaRe = kappa_*Re;
|
||||
|
||||
scalar yp = yPlusLam_;
|
||||
const scalar ryPlusLam = 1.0/yp;
|
||||
|
||||
label iter = 0;
|
||||
scalar yPlusLast = 0.0;
|
||||
|
||||
do
|
||||
{
|
||||
yPlusLast = yp;
|
||||
yp = (kappaRe + yp)/(1.0 + log(E_*yp));
|
||||
|
||||
} while (mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 10);
|
||||
|
||||
yPlus[faceI] = max(0.0, yp);
|
||||
}
|
||||
}
|
||||
|
||||
return tyPlus;
|
||||
}
|
||||
|
||||
|
||||
tmp<scalarField> mutURoughWallFunctionFvPatchScalarField::calcMut() const
|
||||
{
|
||||
const label patchI = patch().index();
|
||||
|
||||
const turbulenceModel& turbModel =
|
||||
db().lookupObject<turbulenceModel>("turbulenceModel");
|
||||
|
||||
const scalarField& y = turbModel.y()[patchI];
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchI];
|
||||
const scalarField& muw = turbModel.mu().boundaryField()[patchI];
|
||||
const fvPatchScalarField& rho = turbModel.rho().boundaryField()[patchI];
|
||||
|
||||
// The flow velocity at the adjacent cell centre
|
||||
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
|
||||
|
||||
tmp<scalarField> tyPlus = calcYPlus(magUp);
|
||||
const scalarField& yPlus = tyPlus();
|
||||
|
||||
tmp<scalarField> tmutw(new scalarField(patch().size(), 0.0));
|
||||
scalarField& mutw = tmutw();
|
||||
|
||||
forAll(yPlus, faceI)
|
||||
{
|
||||
if (yPlus[faceI] > yPlusLam_)
|
||||
{
|
||||
const scalar Re = rho[faceI]*magUp[faceI]*y[faceI]/muw[faceI] + ROOTVSMALL;
|
||||
mutw[faceI] = muw[faceI]*(sqr(yPlus[faceI])/Re - 1);
|
||||
}
|
||||
}
|
||||
|
||||
return tmutw;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
mutURoughWallFunctionFvPatchScalarField::mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const fvPatch& p,
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
)
|
||||
:
|
||||
mutWallFunctionFvPatchScalarField(p, iF),
|
||||
roughnessHeight_(0.0),
|
||||
roughnessConstant_(0.0),
|
||||
roughnessFactor_(0.0)
|
||||
{}
|
||||
|
||||
|
||||
mutURoughWallFunctionFvPatchScalarField::mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const mutURoughWallFunctionFvPatchScalarField& ptf,
|
||||
const fvPatch& p,
|
||||
const DimensionedField<scalar, volMesh>& iF,
|
||||
const fvPatchFieldMapper& mapper
|
||||
)
|
||||
:
|
||||
mutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
|
||||
roughnessHeight_(ptf.roughnessHeight_),
|
||||
roughnessConstant_(ptf.roughnessConstant_),
|
||||
roughnessFactor_(ptf.roughnessFactor_)
|
||||
{}
|
||||
|
||||
|
||||
mutURoughWallFunctionFvPatchScalarField::mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const fvPatch& p,
|
||||
const DimensionedField<scalar, volMesh>& iF,
|
||||
const dictionary& dict
|
||||
)
|
||||
:
|
||||
mutWallFunctionFvPatchScalarField(p, iF, dict),
|
||||
roughnessHeight_(readScalar(dict.lookup("roughnessHeight"))),
|
||||
roughnessConstant_(readScalar(dict.lookup("roughnessConstant"))),
|
||||
roughnessFactor_(readScalar(dict.lookup("roughnessFactor")))
|
||||
{}
|
||||
|
||||
|
||||
mutURoughWallFunctionFvPatchScalarField::mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const mutURoughWallFunctionFvPatchScalarField& rwfpsf
|
||||
)
|
||||
:
|
||||
mutWallFunctionFvPatchScalarField(rwfpsf),
|
||||
roughnessHeight_(rwfpsf.roughnessHeight_),
|
||||
roughnessConstant_(rwfpsf.roughnessConstant_),
|
||||
roughnessFactor_(rwfpsf.roughnessFactor_)
|
||||
{}
|
||||
|
||||
|
||||
mutURoughWallFunctionFvPatchScalarField::mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const mutURoughWallFunctionFvPatchScalarField& rwfpsf,
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
)
|
||||
:
|
||||
mutWallFunctionFvPatchScalarField(rwfpsf, iF),
|
||||
roughnessHeight_(rwfpsf.roughnessHeight_),
|
||||
roughnessConstant_(rwfpsf.roughnessConstant_),
|
||||
roughnessFactor_(rwfpsf.roughnessFactor_)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||
|
||||
tmp<scalarField> mutURoughWallFunctionFvPatchScalarField::yPlus() const
|
||||
{
|
||||
const label patchI = patch().index();
|
||||
|
||||
const turbulenceModel& turbModel =
|
||||
db().lookupObject<turbulenceModel>("turbulenceModel");
|
||||
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchI];
|
||||
const scalarField magUp = mag(Uw.patchInternalField() - Uw);
|
||||
|
||||
return calcYPlus(magUp);
|
||||
}
|
||||
|
||||
|
||||
void mutURoughWallFunctionFvPatchScalarField::write(Ostream& os) const
|
||||
{
|
||||
mutWallFunctionFvPatchScalarField::write(os);
|
||||
writeLocalEntries(os);
|
||||
os.writeKeyword("roughnessHeight")
|
||||
<< roughnessHeight_ << token::END_STATEMENT << nl;
|
||||
os.writeKeyword("roughnessConstant")
|
||||
<< roughnessConstant_ << token::END_STATEMENT << nl;
|
||||
os.writeKeyword("roughnessFactor")
|
||||
<< roughnessFactor_ << token::END_STATEMENT << nl;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
makePatchTypeField
|
||||
(
|
||||
fvPatchScalarField,
|
||||
mutURoughWallFunctionFvPatchScalarField
|
||||
);
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace RASModels
|
||||
} // End namespace compressible
|
||||
} // End namespace Foam
|
||||
|
||||
// ************************************************************************* //
|
|
@ -0,0 +1,237 @@
|
|||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | foam-extend: Open Source CFD
|
||||
\\ / O peration | Version: 4.1
|
||||
\\ / A nd | Web: http://www.foam-extend.org
|
||||
\\/ M anipulation | For copyright notice see file Copyright
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
This file is part of foam-extend.
|
||||
|
||||
foam-extend 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.
|
||||
|
||||
foam-extend 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 foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
Class
|
||||
Foam::compressible::RASModels::mutURoughWallFunctionFvPatchScalarField
|
||||
|
||||
Description
|
||||
This boundary condition provides a turbulent (dynamic) viscosity condition
|
||||
when using wall functions for rough walls, based on velocity.
|
||||
|
||||
Usage
|
||||
\table
|
||||
Property | Description | Required | Default value
|
||||
roughnessHeight | roughness height | yes |
|
||||
roughnessConstant | roughness constanr | yes |
|
||||
roughnessFactor | scaling factor | yes |
|
||||
\endtable
|
||||
|
||||
Example of the boundary condition specification:
|
||||
\verbatim
|
||||
<patchName>
|
||||
{
|
||||
type mutURoughWallFunction;
|
||||
roughnessHeight 1e-5;
|
||||
roughnessConstant 0.5;
|
||||
roughnessFactor 1;
|
||||
}
|
||||
\endverbatim
|
||||
|
||||
See also
|
||||
Foam::mutWallFunctionFvPatchScalarField
|
||||
|
||||
SourceFiles
|
||||
mutURoughWallFunctionFvPatchScalarField.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef mutURoughWallFunctionFvPatchScalarField_H
|
||||
#define mutURoughWallFunctionFvPatchScalarField_H
|
||||
|
||||
#include "mutWallFunctionFvPatchScalarField.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace compressible
|
||||
{
|
||||
namespace RASModels
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class mutURoughWallFunctionFvPatchScalarField Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class mutURoughWallFunctionFvPatchScalarField
|
||||
:
|
||||
public mutWallFunctionFvPatchScalarField
|
||||
{
|
||||
// Private Data
|
||||
|
||||
// Roughness model parameters
|
||||
|
||||
//- Height
|
||||
scalar roughnessHeight_;
|
||||
|
||||
//- Constant
|
||||
scalar roughnessConstant_;
|
||||
|
||||
//- Scale factor
|
||||
scalar roughnessFactor_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Calculate yPLus
|
||||
virtual tmp<scalarField> calcYPlus(const scalarField& magUp) const;
|
||||
|
||||
//- Calculate the turbulence viscosity
|
||||
virtual tmp<scalarField> calcMut() const;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("mutURoughWallFunction");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from patch and internal field
|
||||
mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const fvPatch&,
|
||||
const DimensionedField<scalar, volMesh>&
|
||||
);
|
||||
|
||||
//- Construct from patch, internal field and dictionary
|
||||
mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const fvPatch&,
|
||||
const DimensionedField<scalar, volMesh>&,
|
||||
const dictionary&
|
||||
);
|
||||
|
||||
//- Construct by mapping given
|
||||
// mutURoughWallFunctionFvPatchScalarField
|
||||
// onto a new patch
|
||||
mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const mutURoughWallFunctionFvPatchScalarField&,
|
||||
const fvPatch&,
|
||||
const DimensionedField<scalar, volMesh>&,
|
||||
const fvPatchFieldMapper&
|
||||
);
|
||||
|
||||
//- Construct as copy
|
||||
mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const mutURoughWallFunctionFvPatchScalarField&
|
||||
);
|
||||
|
||||
//- Construct and return a clone
|
||||
virtual tmp<fvPatchScalarField> clone() const
|
||||
{
|
||||
return tmp<fvPatchScalarField>
|
||||
(
|
||||
new mutURoughWallFunctionFvPatchScalarField(*this)
|
||||
);
|
||||
}
|
||||
|
||||
//- Construct as copy setting internal field reference
|
||||
mutURoughWallFunctionFvPatchScalarField
|
||||
(
|
||||
const mutURoughWallFunctionFvPatchScalarField&,
|
||||
const DimensionedField<scalar, volMesh>&
|
||||
);
|
||||
|
||||
//- Construct and return a clone setting internal field reference
|
||||
virtual tmp<fvPatchScalarField> clone
|
||||
(
|
||||
const DimensionedField<scalar, volMesh>& iF
|
||||
) const
|
||||
{
|
||||
return tmp<fvPatchScalarField>
|
||||
(
|
||||
new mutURoughWallFunctionFvPatchScalarField(*this, iF)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
// Member functions
|
||||
|
||||
// Access
|
||||
|
||||
//- Return the roughness height scale
|
||||
const scalar& roughnessHeight() const
|
||||
{
|
||||
return roughnessHeight_;
|
||||
}
|
||||
|
||||
//- Return reference to the roughness height to allow adjustment
|
||||
scalar& roughnessHeight()
|
||||
{
|
||||
return roughnessHeight_;
|
||||
}
|
||||
|
||||
//- Return the roughness constant scale
|
||||
const scalar& roughnessConstant() const
|
||||
{
|
||||
return roughnessConstant_;
|
||||
}
|
||||
|
||||
//- Return reference to the roughness constant to allow adjustment
|
||||
scalar& roughnessConstant()
|
||||
{
|
||||
return roughnessConstant_;
|
||||
}
|
||||
|
||||
//- Return the roughness scale factor
|
||||
const scalar& roughnessFudgeFactor() const
|
||||
{
|
||||
return roughnessFactor_;
|
||||
}
|
||||
|
||||
//- Return reference to the roughness scale factor to allow
|
||||
// adjustment
|
||||
scalar& roughnessFactor()
|
||||
{
|
||||
return roughnessFactor_;
|
||||
}
|
||||
|
||||
|
||||
// Evaluation functions
|
||||
|
||||
//- Calculate and return the yPlus at the boundary
|
||||
virtual tmp<scalarField> yPlus() const;
|
||||
|
||||
|
||||
// I-O
|
||||
|
||||
//- Write
|
||||
virtual void write(Ostream& os) const;
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace RASModels
|
||||
} // End namespace compressible
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
Reference in a new issue