epsilonmWallFunction: New wall-function specifically for the mixtureKEpsilon model

epsilonm is obtained by combining epsilon.gas and epsilon.liquid in a two-phase
system, each of which will apply the epsilonWallFunction at walls; the
epsilonmWallFunction propagates the resulting wall epsilonm into the near-wall
cells.

If the 0/epsilonm file is provided the epsilonmWallFunction should be specified
for walls, if the 0/epsilonm file is not provided it will be generated
automatically and the epsilonmWallFunction applied to walls for which the
epsilonWallFunction is specified in the epsilon.liquid file.
This commit is contained in:
Henry Weller
2022-04-20 18:48:35 +01:00
parent a4a03db574
commit 3bac211785
6 changed files with 331 additions and 3 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -204,6 +204,7 @@ public:
);
}
//- Destructor
virtual ~epsilonWallFunctionFvPatchScalarField()
{}

View File

@ -1,4 +1,6 @@
phaseCompressibleMomentumTransportModel.C
phaseCompressibleMomentumTransportModels.C
derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonmWallFunction/epsilonmWallFunctionFvPatchScalarField.C
LIB = $(FOAM_LIBBIN)/libphaseCompressibleMomentumTransportModels

View File

@ -31,6 +31,8 @@ License
#include "dispersedVirtualMassModel.H"
#include "fixedValueFvPatchFields.H"
#include "inletOutletFvPatchFields.H"
#include "epsilonWallFunctionFvPatchScalarField.H"
#include "epsilonmWallFunctionFvPatchScalarField.H"
#include "fvmSup.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -187,7 +189,11 @@ wordList mixtureKEpsilon<BasicMomentumTransportModel>::epsilonBoundaryTypes
forAll(ebf, patchi)
{
if (isA<fixedValueFvPatchScalarField>(ebf[patchi]))
if (isType<epsilonWallFunctionFvPatchScalarField>(ebf[patchi]))
{
ebt[patchi] = epsilonmWallFunctionFvPatchScalarField::typeName;
}
else if (isA<fixedValueFvPatchScalarField>(ebf[patchi]))
{
ebt[patchi] = fixedValueFvPatchScalarField::typeName;
}

View File

@ -0,0 +1,170 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\/ 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 "epsilonmWallFunctionFvPatchScalarField.H"
#include "fvMatrix.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
Foam::scalar Foam::epsilonmWallFunctionFvPatchScalarField::tolerance_ = 1e-5;
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::epsilonmWallFunctionFvPatchScalarField::
epsilonmWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchField<scalar>(p, iF)
{}
Foam::epsilonmWallFunctionFvPatchScalarField::
epsilonmWallFunctionFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchField<scalar>(p, iF, dict)
{
// Apply zero-gradient condition on start-up
this->operator==(patchInternalField());
}
Foam::epsilonmWallFunctionFvPatchScalarField::
epsilonmWallFunctionFvPatchScalarField
(
const epsilonmWallFunctionFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField<scalar>(ptf, p, iF, mapper)
{}
Foam::epsilonmWallFunctionFvPatchScalarField::
epsilonmWallFunctionFvPatchScalarField
(
const epsilonmWallFunctionFvPatchScalarField& ewfpsf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchField<scalar>(ewfpsf, iF)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::epsilonmWallFunctionFvPatchScalarField::manipulateMatrix
(
fvMatrix<scalar>& matrix
)
{
if (manipulatedMatrix())
{
return;
}
matrix.setValues(patch().faceCells(), patchInternalField()());
fvPatchField<scalar>::manipulateMatrix(matrix);
}
void Foam::epsilonmWallFunctionFvPatchScalarField::manipulateMatrix
(
fvMatrix<scalar>& matrix,
const Field<scalar>& weights
)
{
if (manipulatedMatrix())
{
return;
}
DynamicList<label> constraintCells(weights.size());
DynamicList<scalar> constraintEpsilon(weights.size());
const labelUList& faceCells = patch().faceCells();
const DimensionedField<scalar, volMesh>& epsilon
= internalField();
label nConstrainedCells = 0;
forAll(weights, facei)
{
// Only set the values if the weights are > tolerance
if (weights[facei] > tolerance_)
{
nConstrainedCells++;
label celli = faceCells[facei];
constraintCells.append(celli);
constraintEpsilon.append(epsilon[celli]);
}
}
if (debug)
{
Pout<< "Patch: " << patch().name()
<< ": number of constrained cells = " << nConstrainedCells
<< " out of " << patch().size()
<< endl;
}
matrix.setValues
(
constraintCells,
scalarField(constraintEpsilon)
);
fvPatchField<scalar>::manipulateMatrix(matrix);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
epsilonmWallFunctionFvPatchScalarField
);
}
// ************************************************************************* //

View File

@ -0,0 +1,149 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\/ 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/>.
Class
Foam::epsilonmWallFunctionFvPatchScalarField
Description
This boundary condition provides a turbulence dissipation wall constraint
for the Foam::mixtureKEpsilon model.
SourceFiles
epsilonmWallFunctionFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef epsilonmWallFunctionFvPatchScalarField_H
#define epsilonmWallFunctionFvPatchScalarField_H
#include "fixedValueFvPatchField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class epsilonmWallFunctionFvPatchScalarField Declaration
\*---------------------------------------------------------------------------*/
class epsilonmWallFunctionFvPatchScalarField
:
public fixedValueFvPatchField<scalar>
{
// Private member data
//- Tolerance used in weighted calculations
static scalar tolerance_;
public:
//- Runtime type information
TypeName("epsilonmWallFunction");
// Constructors
//- Construct from patch and internal field
epsilonmWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
epsilonmWallFunctionFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given
// epsilonmWallFunctionFvPatchScalarField
// onto a new patch
epsilonmWallFunctionFvPatchScalarField
(
const epsilonmWallFunctionFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Disallow copy without setting internal field reference
epsilonmWallFunctionFvPatchScalarField
(
const epsilonmWallFunctionFvPatchScalarField&
) = delete;
//- Copy constructor setting internal field reference
epsilonmWallFunctionFvPatchScalarField
(
const epsilonmWallFunctionFvPatchScalarField&,
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 epsilonmWallFunctionFvPatchScalarField(*this, iF)
);
}
//- Destructor
virtual ~epsilonmWallFunctionFvPatchScalarField()
{}
// Member Functions
// Evaluation functions
//- Manipulate matrix
virtual void manipulateMatrix(fvMatrix<scalar>& matrix);
//- Manipulate matrix with given weights
virtual void manipulateMatrix
(
fvMatrix<scalar>& matrix,
const scalarField& weights
);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -35,7 +35,7 @@ boundaryField
walls
{
type zeroGradient;
type epsilonmWallFunction;
value $internalField;
}