The patch-specific mapper interfaces, fvPatchFieldMapper and pointPatchFieldMapper, have been removed as they did not do anything. Patch mapping constructors and functions now take a basic fieldMapper reference. An fvPatchFieldMapper.H header has been provided to aid backwards compatability so that existing custom boundary conditions continue to compile.
189 lines
5.0 KiB
C++
189 lines
5.0 KiB
C++
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration | Website: https://openfoam.org
|
|
\\ / A nd | Copyright (C) 2014-2023 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 "JohnsonJacksonParticleSlipFvPatchVectorField.H"
|
|
#include "addToRunTimeSelectionTable.H"
|
|
#include "kineticTheoryModel.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
makePatchTypeField
|
|
(
|
|
fvPatchVectorField,
|
|
JohnsonJacksonParticleSlipFvPatchVectorField
|
|
);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
|
|
JohnsonJacksonParticleSlipFvPatchVectorField
|
|
(
|
|
const fvPatch& p,
|
|
const DimensionedField<vector, volMesh>& iF,
|
|
const dictionary& dict
|
|
)
|
|
:
|
|
partialSlipFvPatchVectorField(p, iF),
|
|
specularityCoefficient_
|
|
(
|
|
"specularityCoefficient",
|
|
dimless,
|
|
dict.lookup("specularityCoefficient")
|
|
)
|
|
{
|
|
if
|
|
(
|
|
(specularityCoefficient_.value() < 0)
|
|
|| (specularityCoefficient_.value() > 1)
|
|
)
|
|
{
|
|
FatalErrorInFunction
|
|
<< "The specularity coefficient has to be between 0 and 1"
|
|
<< abort(FatalError);
|
|
}
|
|
|
|
fvPatchVectorField::operator=
|
|
(
|
|
vectorField("value", dict, p.size())
|
|
);
|
|
}
|
|
|
|
|
|
Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
|
|
JohnsonJacksonParticleSlipFvPatchVectorField
|
|
(
|
|
const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
|
|
const fvPatch& p,
|
|
const DimensionedField<vector, volMesh>& iF,
|
|
const fieldMapper& mapper
|
|
)
|
|
:
|
|
partialSlipFvPatchVectorField(ptf, p, iF, mapper),
|
|
specularityCoefficient_(ptf.specularityCoefficient_)
|
|
{}
|
|
|
|
|
|
Foam::JohnsonJacksonParticleSlipFvPatchVectorField::
|
|
JohnsonJacksonParticleSlipFvPatchVectorField
|
|
(
|
|
const JohnsonJacksonParticleSlipFvPatchVectorField& ptf,
|
|
const DimensionedField<vector, volMesh>& iF
|
|
)
|
|
:
|
|
partialSlipFvPatchVectorField(ptf, iF),
|
|
specularityCoefficient_(ptf.specularityCoefficient_)
|
|
{}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
|
|
|
void Foam::JohnsonJacksonParticleSlipFvPatchVectorField::updateCoeffs()
|
|
{
|
|
if (updated())
|
|
{
|
|
return;
|
|
}
|
|
|
|
// lookup the fluid model and the phase
|
|
const phaseSystem& fluid =
|
|
db().lookupObject<phaseSystem>(phaseSystem::propertiesName);
|
|
|
|
const phaseModel& phase
|
|
(
|
|
fluid.phases()[internalField().group()]
|
|
);
|
|
|
|
// lookup all the fields on this patch
|
|
const fvPatchScalarField& alpha
|
|
(
|
|
patch().lookupPatchField<volScalarField, scalar>
|
|
(
|
|
phase.volScalarField::name()
|
|
)
|
|
);
|
|
|
|
const fvPatchScalarField& gs0
|
|
(
|
|
patch().lookupPatchField<volScalarField, scalar>
|
|
(
|
|
IOobject::groupName
|
|
(
|
|
Foam::typedName<RASModels::kineticTheoryModel>("gs0"),
|
|
phase.name()
|
|
)
|
|
)
|
|
);
|
|
|
|
const scalarField nu
|
|
(
|
|
patch().lookupPatchField<volScalarField, scalar>
|
|
(
|
|
IOobject::groupName("nut", phase.name())
|
|
)
|
|
);
|
|
|
|
word ThetaName(IOobject::groupName("Theta", phase.name()));
|
|
|
|
const fvPatchScalarField& Theta
|
|
(
|
|
db().foundObject<volScalarField>(ThetaName)
|
|
? patch().lookupPatchField<volScalarField, scalar>(ThetaName)
|
|
: alpha
|
|
);
|
|
|
|
// calculate the slip value fraction
|
|
scalarField c
|
|
(
|
|
constant::mathematical::pi
|
|
*alpha
|
|
*gs0
|
|
*specularityCoefficient_.value()
|
|
*sqrt(3*Theta)
|
|
/max(6*nu*phase.alphaMax(), small)
|
|
);
|
|
|
|
this->valueFraction() = c/(c + patch().deltaCoeffs());
|
|
|
|
partialSlipFvPatchVectorField::updateCoeffs();
|
|
}
|
|
|
|
|
|
void Foam::JohnsonJacksonParticleSlipFvPatchVectorField::write
|
|
(
|
|
Ostream& os
|
|
) const
|
|
{
|
|
fvPatchVectorField::write(os);
|
|
writeEntry(os, "specularityCoefficient", specularityCoefficient_);
|
|
writeEntry(os, "value", *this);
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|