Refactored nutkRoughWallFunction according to Vanilla

Derives from nutkWallFunctionFvPatchScalarField
This commit is contained in:
Vuko Vukcevic 2018-06-19 11:15:00 +02:00
parent c4db345abe
commit 0c17c8057e
2 changed files with 466 additions and 0 deletions

View file

@ -0,0 +1,249 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "nutkRoughWallFunctionFvPatchScalarField.H"
#include "RASModel.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
scalar nutkRoughWallFunctionFvPatchScalarField::fnRough
(
const scalar KsPlus,
const scalar Cs
) const
{
// Return fn based on non-dimensional roughness height
if (KsPlus < 90.0)
{
return pow
(
(KsPlus - 2.25)/87.75 + Cs*KsPlus,
sin(0.4258*(log(KsPlus) - 0.811))
);
}
else
{
return (1.0 + Cs*KsPlus);
}
}
tmp<scalarField> nutkRoughWallFunctionFvPatchScalarField::calcNut() const
{
const label patchI = patch().index();
const turbulenceModel& turbModel =
db().lookupObject<turbulenceModel>("turbulenceModel");
const scalarField& y = turbModel.y()[patchI];
const tmp<volScalarField> tk = turbModel.k();
const volScalarField& k = tk();
const scalarField& nuw = turbModel.nu().boundaryField()[patchI];
const scalar Cmu25 = pow025(Cmu_);
tmp<scalarField> tnutw(new scalarField(*this));
scalarField& nutw = tnutw();
const unallocLabelList& fc = patch().faceCells();
forAll(nutw, faceI)
{
const label faceCellI = fc[faceI];
const scalar uStar = Cmu25*sqrt(k[faceCellI]);
const scalar yPlus = uStar*y[faceI]/nuw[faceI];
const scalar KsPlus = uStar*Ks_[faceI]/nuw[faceI];
scalar Edash = E_;
if (KsPlus > 2.25)
{
Edash /= fnRough(KsPlus, Cs_[faceI]);
}
if (yPlus > yPlusLam_)
{
scalar limitingNutw = max(nutw[faceI], nuw[faceI]);
// To avoid oscillations limit the change in the wall viscosity
// which is particularly important if it temporarily becomes zero
nutw[faceI] =
max
(
min
(
nuw[faceI]
*(yPlus*kappa_/log(max(Edash*yPlus, 1+1e-4)) - 1),
2*limitingNutw
), 0.5*limitingNutw
);
}
if (debug)
{
Info<< "yPlus = " << yPlus
<< ", KsPlus = " << KsPlus
<< ", Edash = " << Edash
<< ", nutw = " << nutw[faceI]
<< endl;
}
}
return tnutw;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
nutkWallFunctionFvPatchScalarField(p, iF),
Ks_(p.size(), 0.0),
Cs_(p.size(), 0.0)
{}
nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
(
const nutkRoughWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
Ks_(ptf.Ks_, mapper),
Cs_(ptf.Cs_, mapper)
{}
nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
nutkWallFunctionFvPatchScalarField(p, iF, dict),
Ks_("Ks", dict, p.size()),
Cs_("Cs", dict, p.size())
{}
nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
(
const nutkRoughWallFunctionFvPatchScalarField& rwfpsf
)
:
nutkWallFunctionFvPatchScalarField(rwfpsf),
Ks_(rwfpsf.Ks_),
Cs_(rwfpsf.Cs_)
{}
nutkRoughWallFunctionFvPatchScalarField::nutkRoughWallFunctionFvPatchScalarField
(
const nutkRoughWallFunctionFvPatchScalarField& rwfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
nutkWallFunctionFvPatchScalarField(rwfpsf, iF),
Ks_(rwfpsf.Ks_),
Cs_(rwfpsf.Cs_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void nutkRoughWallFunctionFvPatchScalarField::autoMap
(
const fvPatchFieldMapper& m
)
{
nutkWallFunctionFvPatchScalarField::autoMap(m);
Ks_.autoMap(m);
Cs_.autoMap(m);
}
void nutkRoughWallFunctionFvPatchScalarField::rmap
(
const fvPatchScalarField& ptf,
const labelList& addr
)
{
nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
const nutkRoughWallFunctionFvPatchScalarField& nrwfpsf =
refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
Ks_.rmap(nrwfpsf.Ks_, addr);
Cs_.rmap(nrwfpsf.Cs_, addr);
}
void nutkRoughWallFunctionFvPatchScalarField::write(Ostream& os) const
{
fvPatchField<scalar>::write(os);
writeLocalEntries(os);
Cs_.writeEntry("Cs", os);
Ks_.writeEntry("Ks", os);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
makePatchTypeField
(
fvPatchScalarField,
nutkRoughWallFunctionFvPatchScalarField
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// ************************************************************************* //

View file

@ -0,0 +1,217 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::incompressible::RASModels::nutkRoughWallFunctionFvPatchScalarField
Description
Boundary condition for turbulent (kinematic) viscosity when using wall
functions for rough walls.
Manipulates the E parameter to account for roughness effects, based on
KsPlus.
- roughness height = sand-grain roughness (0 for smooth walls)
- roughness constant = 0.5-1.0 (0.5 default)
Usage
\table
Property | Description | Required | Default value
Ks | sand-grain roughness height | yes |
Cs | roughness constant | yes |
\endtable
Example of the boundary condition specification:
\verbatim
<patchName>
{
type nutkRoughWallFunction;
Ks uniform 0;
Cs uniform 0.5;
}
\endverbatim
SourceFiles
nutkRoughWallFunctionFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef nutkRoughWallFunctionFvPatchScalarField_H
#define nutkRoughWallFunctionFvPatchScalarField_H
#include "nutkWallFunctionFvPatchScalarField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace incompressible
{
namespace RASModels
{
/*---------------------------------------------------------------------------*\
Class nutkRoughWallFunctionFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
class nutkRoughWallFunctionFvPatchScalarField
:
public nutkWallFunctionFvPatchScalarField
{
protected:
// Protected data
//- Roughness height
scalarField Ks_;
//- Roughness constant
scalarField Cs_;
// Protected Member Functions
//- Compute the roughness function
virtual scalar fnRough(const scalar KsPlus, const scalar Cs) const;
//- Calculate the turbulence viscosity
virtual tmp<scalarField> calcNut() const;
public:
//- Runtime type information
TypeName("nutkRoughWallFunction");
// Constructors
//- Construct from patch and internal field
nutkRoughWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
nutkRoughWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given
// nutkRoughWallFunctionFvPatchScalarField
// onto a new patch
nutkRoughWallFunctionFvPatchScalarField
(
const nutkRoughWallFunctionFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
nutkRoughWallFunctionFvPatchScalarField
(
const nutkRoughWallFunctionFvPatchScalarField&
);
//- Construct and return a clone
virtual tmp<fvPatchScalarField> clone() const
{
return tmp<fvPatchScalarField>
(
new nutkRoughWallFunctionFvPatchScalarField(*this)
);
}
//- Construct as copy setting internal field reference
nutkRoughWallFunctionFvPatchScalarField
(
const nutkRoughWallFunctionFvPatchScalarField&,
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 nutkRoughWallFunctionFvPatchScalarField(*this, iF)
);
}
// Member functions
// Acces functions
// Return Ks
scalarField& Ks()
{
return Ks_;
}
// Return Cs
scalarField& Cs()
{
return Cs_;
}
// Mapping functions
//- Map (and resize as needed) from self given a mapping object
virtual void autoMap(const fvPatchFieldMapper&);
//- Reverse map the given fvPatchField onto this fvPatchField
virtual void rmap
(
const fvPatchScalarField&,
const labelList&
);
// I-O
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace RASModels
} // End namespace incompressible
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //