mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: nutWallFunctions - provided option to specify alternative velocity field
This commit is contained in:
@ -114,7 +114,7 @@ tmp<scalarField> nutLowReWallFunctionFvPatchScalarField::yPlus() const
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
|
||||
return y*sqrt(nuw*mag(Uw.snGrad()))/nuw;
|
||||
}
|
||||
|
||||
@ -44,7 +44,7 @@ Foam::nutUBlendedWallFunctionFvPatchScalarField::calcNut() const
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magGradU(mag(Uw.snGrad()));
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
@ -80,7 +80,7 @@ Foam::nutUBlendedWallFunctionFvPatchScalarField::calcUTau
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
const vectorField n(patch().nf());
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
vectorField Up(Uw.patchInternalField() - Uw);
|
||||
Up -= n*(n & Up);
|
||||
const scalarField magUp(mag(Up));
|
||||
@ -197,7 +197,7 @@ Foam::nutUBlendedWallFunctionFvPatchScalarField::yPlus() const
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magGradU(mag(Uw.snGrad()));
|
||||
|
||||
return y*calcUTau(magGradU)/nuw;
|
||||
|
||||
@ -51,7 +51,7 @@ tmp<scalarField> nutURoughWallFunctionFvPatchScalarField::calcNut() const
|
||||
)
|
||||
);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
@ -314,7 +314,7 @@ tmp<scalarField> nutURoughWallFunctionFvPatchScalarField::yPlus() const
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
tmp<scalarField> magUp = mag(Uw.patchInternalField() - Uw);
|
||||
|
||||
return calcYPlus(magUp());
|
||||
|
||||
@ -50,7 +50,7 @@ tmp<scalarField> nutUSpaldingWallFunctionFvPatchScalarField::calcNut() const
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magGradU(mag(Uw.snGrad()));
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
@ -120,7 +120,7 @@ tmp<scalarField> nutUSpaldingWallFunctionFvPatchScalarField::calcUTau
|
||||
);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
|
||||
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
@ -327,7 +327,7 @@ tmp<scalarField> nutUSpaldingWallFunctionFvPatchScalarField::yPlus() const
|
||||
)
|
||||
);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
|
||||
@ -51,7 +51,7 @@ tmp<scalarField> nutUTabulatedWallFunctionFvPatchScalarField::calcNut() const
|
||||
)
|
||||
);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
|
||||
const scalarField magGradU(mag(Uw.snGrad()));
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
@ -193,7 +193,7 @@ tmp<scalarField> nutUTabulatedWallFunctionFvPatchScalarField::yPlus() const
|
||||
)
|
||||
);
|
||||
const scalarField& y = turbModel.y()[patchi];
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
|
||||
@ -50,7 +50,7 @@ tmp<scalarField> nutUWallFunctionFvPatchScalarField::calcNut() const
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
|
||||
const tmp<scalarField> tnuw = turbModel.nu(patchi);
|
||||
const scalarField& nuw = tnuw();
|
||||
@ -187,7 +187,7 @@ tmp<scalarField> nutUWallFunctionFvPatchScalarField::yPlus() const
|
||||
internalField().group()
|
||||
)
|
||||
);
|
||||
const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
|
||||
const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
|
||||
const scalarField magUp(mag(Uw.patchInternalField() - Uw));
|
||||
|
||||
return calcYPlus(magUp);
|
||||
|
||||
@ -29,6 +29,7 @@ License
|
||||
#include "fvPatchFieldMapper.H"
|
||||
#include "volFields.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "turbulenceModel.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
@ -54,11 +55,32 @@ void Foam::nutWallFunctionFvPatchScalarField::checkType()
|
||||
}
|
||||
|
||||
|
||||
const Foam::volVectorField& Foam::nutWallFunctionFvPatchScalarField::U
|
||||
(
|
||||
const turbulenceModel& turb
|
||||
) const
|
||||
{
|
||||
if (UName_ == word::null)
|
||||
{
|
||||
return turb.U();
|
||||
}
|
||||
else
|
||||
{
|
||||
return db().lookupObject<volVectorField>(UName_);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void Foam::nutWallFunctionFvPatchScalarField::writeLocalEntries
|
||||
(
|
||||
Ostream& os
|
||||
) const
|
||||
{
|
||||
if (UName_ != word::null)
|
||||
{
|
||||
os.writeEntry("U", UName_);
|
||||
}
|
||||
|
||||
os.writeEntry("Cmu", Cmu_);
|
||||
os.writeEntry("kappa", kappa_);
|
||||
os.writeEntry("E", E_);
|
||||
@ -74,6 +96,7 @@ Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(p, iF),
|
||||
UName_(word::null),
|
||||
Cmu_(0.09),
|
||||
kappa_(0.41),
|
||||
E_(9.8),
|
||||
@ -92,6 +115,7 @@ Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
|
||||
UName_(ptf.UName_),
|
||||
Cmu_(ptf.Cmu_),
|
||||
kappa_(ptf.kappa_),
|
||||
E_(ptf.E_),
|
||||
@ -109,6 +133,7 @@ Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(p, iF, dict),
|
||||
UName_(dict.lookupOrDefault<word>("U", word::null)),
|
||||
Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
|
||||
kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
|
||||
E_(dict.lookupOrDefault<scalar>("E", 9.8)),
|
||||
@ -124,6 +149,7 @@ Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(wfpsf),
|
||||
UName_(wfpsf.UName_),
|
||||
Cmu_(wfpsf.Cmu_),
|
||||
kappa_(wfpsf.kappa_),
|
||||
E_(wfpsf.E_),
|
||||
@ -140,6 +166,7 @@ Foam::nutWallFunctionFvPatchScalarField::nutWallFunctionFvPatchScalarField
|
||||
)
|
||||
:
|
||||
fixedValueFvPatchScalarField(wfpsf, iF),
|
||||
UName_(wfpsf.UName_),
|
||||
Cmu_(wfpsf.Cmu_),
|
||||
kappa_(wfpsf.kappa_),
|
||||
E_(wfpsf.E_),
|
||||
|
||||
@ -76,6 +76,8 @@ SourceFiles
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
class turbulenceModel;
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class nutWallFunctionFvPatchScalarField Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -88,6 +90,11 @@ protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Name of velocity field
|
||||
// Defult is null (not specified) in which case the velocity is
|
||||
// retrieved from the turbulence model
|
||||
word UName_;
|
||||
|
||||
//- Cmu coefficient
|
||||
scalar Cmu_;
|
||||
|
||||
@ -103,6 +110,10 @@ protected:
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Helper to return the velocity field either from the turbulence
|
||||
//- model (default) or the mesh database
|
||||
virtual const volVectorField& U(const turbulenceModel& turb) const;
|
||||
|
||||
//- Check the type of the patch
|
||||
virtual void checkType();
|
||||
|
||||
|
||||
Reference in New Issue
Block a user