Files
OpenFOAM-12/src/finiteVolume/fields/fvPatchFields/derived/flowRateInletVelocity/flowRateInletVelocityFvPatchVectorField.C
Henry Weller 8bb48df87f flowRateInletVelocityFvPatchVectorField: Added optional profile entry to specify the velocity profile
The unreliable extrapolateProfile option has been replaced by the more flexible
and reliable profile option which allows the velocity profile to be specified as
a Function1 of the normalised distance to the wall.  To simplify the
specification of the most common velocity profiles the new laminarBL (quadratic
profile) and turbulentBL (1/7th power law) Function1s are provided.

In addition to the new profile option the flow rate can now be specified as a
meanVelocity, volumetricFlowRate or massFlowRate, all of which are Function1s of
time.

The following tutorials have been updated to use the laminarBL profile:
    multiphase/multiphaseEulerFoam/laminar/titaniaSynthesis
    multiphase/multiphaseEulerFoam/laminar/titaniaSynthesisSurface

The following tutorials have been updated to use the turbulentBL profile:
    combustion/reactingFoam/Lagrangian/verticalChannel
    combustion/reactingFoam/Lagrangian/verticalChannelLTS
    combustion/reactingFoam/Lagrangian/verticalChannelSteady
    compressible/rhoPimpleFoam/RAS/angledDuct
    compressible/rhoPimpleFoam/RAS/angledDuctLTS
    compressible/rhoPimpleFoam/RAS/squareBendLiq
    compressible/rhoPorousSimpleFoam/angledDuctImplicit
    compressible/rhoSimpleFoam/angledDuctExplicitFixedCoeff
    compressible/rhoSimpleFoam/squareBend
    compressible/rhoSimpleFoam/squareBendLiq
    heatTransfer/chtMultiRegionFoam/shellAndTubeHeatExchanger
    heatTransfer/chtMultiRegionFoam/shellAndTubeHeatExchanger
    incompressible/porousSimpleFoam/angledDuctImplicit
    incompressible/porousSimpleFoam/straightDuctImplicit
    multiphase/interFoam/RAS/angledDuct

Class
    Foam::flowRateInletVelocityFvPatchVectorField

Description
    Velocity inlet boundary condition creating a velocity field with
    optionally specified profile normal to the patch adjusted to match the
    specified mass flow rate, volumetric flow rate or mean velocity.

    For a mass-based flux:
    - the flow rate should be provided in kg/s
    - if \c rho is "none" the flow rate is in m3/s
    - otherwise \c rho should correspond to the name of the density field
    - if the density field cannot be found in the database, the user must
      specify the inlet density using the \c rhoInlet entry

    For a volumetric-based flux:
    - the flow rate is in m3/s

Usage
    \table
        Property     | Description             | Required    | Default value
        massFlowRate | Mass flow rate [kg/s]   | no          |
        volumetricFlowRate | Volumetric flow rate [m^3/s]| no |
        meanVelocity | Mean velocity [m/s]| no |
        profile      | Velocity profile        | no          |
        rho          | Density field name      | no          | rho
        rhoInlet     | Inlet density           | no          |
        alpha        | Volume fraction field name | no       |
    \endtable

    Example of the boundary condition specification for a volumetric flow rate:
    \verbatim
    <patchName>
    {
        type                flowRateInletVelocity;
        volumetricFlowRate  0.2;
        profile             laminarBL;
    }
    \endverbatim

    Example of the boundary condition specification for a mass flow rate:
     \verbatim
    <patchName>
    {
        type                flowRateInletVelocity;
        massFlowRate        0.2;
        profile             turbulentBL;
        rho                 rho;
        rhoInlet            1.0;
    }
    \endverbatim

    Example of the boundary condition specification for a volumetric flow rate:
    \verbatim
    <patchName>
    {
        type                flowRateInletVelocity;
        meanVelocity        5;
        profile             turbulentBL;
    }
    \endverbatim

    The \c volumetricFlowRate, \c massFlowRate or \c meanVelocity entries are
    \c Function1 of time, see Foam::Function1s.

    The \c profile entry is a \c Function1 of the normalised distance to the
    wall.  Any suitable Foam::Function1s can be used including
    Foam::Function1s::codedFunction1 but Foam::Function1s::laminarBL and
    Foam::Function1s::turbulentBL have been created specifically for this
    purpose and are likely to be appropriate for most cases.

Note
    - \c rhoInlet is required for the case of a mass flow rate, where the
      density field is not available at start-up
    - The value is positive into the domain (as an inlet)
    - May not work correctly for transonic inlets
    - Strange behaviour with potentialFoam since the U equation is not solved

See also
    Foam::fixedValueFvPatchField
    Foam::Function1s::laminarBL
    Foam::Function1s::turbulentBL
    Foam::Function1s
    Foam::flowRateOutletVelocityFvPatchVectorField
2022-01-24 19:10:39 +00:00

335 lines
8.3 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-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 "flowRateInletVelocityFvPatchVectorField.H"
#include "addToRunTimeSelectionTable.H"
#include "volFields.H"
#include "one.H"
#include "patchPatchDist.H"
#include "wallPolyPatch.H"
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
void Foam::flowRateInletVelocityFvPatchVectorField::setWallDist()
{
if (profile_.valid())
{
const labelHashSet otherPatchIDs
(
patch().patch().boundaryMesh().findPatchIDs<wallPolyPatch>()
);
const patchPatchDist pwd(patch().patch(), otherPatchIDs);
y_ = pwd/gMax(pwd);
}
area_ = gSum(patch().magSf());
}
template<class ScaleType, class AlphaType, class RhoType>
void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
(
const ScaleType& scale,
const AlphaType& alpha,
const RhoType& rho
)
{
const scalarField profile(this->profile());
const scalar avgU =
-(scale*flowRate_->value(db().time().userTimeValue()))
/gSum(alpha*rho*profile*patch().magSf());
operator==(avgU*profile*patch().nf());
}
template<class AlphaType>
void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
(
const AlphaType& alpha
)
{
if (meanVelocity_)
{
updateValues(area_, alpha, one());
}
else if (volumetric_ || rhoName_ == "none")
{
updateValues(one(), alpha, one());
}
else
{
// Mass flow-rate
if (db().foundObject<volScalarField>(rhoName_))
{
const fvPatchField<scalar>& rhop =
patch().lookupPatchField<volScalarField, scalar>(rhoName_);
updateValues(one(), alpha, rhop);
}
else
{
// Use constant density
if (rhoInlet_ < 0)
{
FatalErrorInFunction
<< "Did not find registered density field " << rhoName_
<< " and no constant density 'rhoInlet' specified"
<< exit(FatalError);
}
updateValues(one(), alpha, rhoInlet_);
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::flowRateInletVelocityFvPatchVectorField::
flowRateInletVelocityFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF
)
:
fixedValueFvPatchField<vector>(p, iF),
flowRate_(),
profile_(),
meanVelocity_(false),
volumetric_(false),
rhoName_("rho"),
rhoInlet_(0),
alphaName_(word::null),
area_(0)
{}
Foam::flowRateInletVelocityFvPatchVectorField::
flowRateInletVelocityFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchField<vector>(p, iF, dict, false),
rhoInlet_(dict.lookupOrDefault<scalar>("rhoInlet", -vGreat)),
alphaName_(dict.lookupOrDefault<word>("alpha", word::null))
{
if (dict.found("meanVelocity"))
{
meanVelocity_ = true;
volumetric_ = false;
flowRate_ = Function1<scalar>::New("meanVelocity", dict);
}
else if (dict.found("volumetricFlowRate"))
{
meanVelocity_ = false;
volumetric_ = true;
flowRate_ = Function1<scalar>::New("volumetricFlowRate", dict);
}
else if (dict.found("massFlowRate"))
{
meanVelocity_ = false;
volumetric_ = false;
flowRate_ = Function1<scalar>::New("massFlowRate", dict);
rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
}
else
{
FatalIOErrorInFunction
(
dict
) << "Please supply 'meanVelocity', 'volumetricFlowRate' or"
<< " 'massFlowRate' and 'rho'" << exit(FatalIOError);
}
if (dict.found("profile"))
{
profile_ = Function1<scalar>::New("profile", dict);
}
setWallDist();
// Value field require if mass based
if (dict.found("value"))
{
fvPatchField<vector>::operator=
(
vectorField("value", dict, p.size())
);
}
else
{
evaluate(Pstream::commsTypes::blocking);
}
}
Foam::flowRateInletVelocityFvPatchVectorField::
flowRateInletVelocityFvPatchVectorField
(
const flowRateInletVelocityFvPatchVectorField& ptf,
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
flowRate_(ptf.flowRate_, false),
profile_(ptf.profile_, false),
meanVelocity_(ptf.meanVelocity_),
volumetric_(ptf.volumetric_),
rhoName_(ptf.rhoName_),
rhoInlet_(ptf.rhoInlet_),
alphaName_(ptf.alphaName_)
{
setWallDist();
}
Foam::flowRateInletVelocityFvPatchVectorField::
flowRateInletVelocityFvPatchVectorField
(
const flowRateInletVelocityFvPatchVectorField& ptf,
const DimensionedField<vector, volMesh>& iF
)
:
fixedValueFvPatchField<vector>(ptf, iF),
flowRate_(ptf.flowRate_, false),
profile_(ptf.profile_, false),
meanVelocity_(ptf.meanVelocity_),
volumetric_(ptf.volumetric_),
rhoName_(ptf.rhoName_),
rhoInlet_(ptf.rhoInlet_),
alphaName_(ptf.alphaName_),
y_(ptf.y_),
area_(ptf.area_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::flowRateInletVelocityFvPatchVectorField::autoMap
(
const fvPatchFieldMapper& m
)
{
fixedValueFvPatchVectorField::autoMap(m);
setWallDist();
}
void Foam::flowRateInletVelocityFvPatchVectorField::rmap
(
const fvPatchVectorField& ptf,
const labelList& addr
)
{
fixedValueFvPatchVectorField::rmap(ptf, addr);
const flowRateInletVelocityFvPatchVectorField& tiptf =
refCast<const flowRateInletVelocityFvPatchVectorField>(ptf);
if (profile_.valid())
{
y_.rmap(tiptf.y_, addr);
}
}
Foam::tmp<Foam::scalarField>
Foam::flowRateInletVelocityFvPatchVectorField::profile()
{
if (profile_.valid())
{
return profile_->value(y_);
}
else
{
return tmp<scalarField>(new scalarField(size(), scalar(1)));
}
}
void Foam::flowRateInletVelocityFvPatchVectorField::updateCoeffs()
{
if (updated())
{
return;
}
if (alphaName_ != word::null)
{
const fvPatchField<scalar>& alphap =
patch().lookupPatchField<volScalarField, scalar>(alphaName_);
updateValues(alphap);
}
else
{
updateValues(one());
}
fixedValueFvPatchVectorField::updateCoeffs();
}
void Foam::flowRateInletVelocityFvPatchVectorField::write(Ostream& os) const
{
fvPatchField<vector>::write(os);
writeEntry(os, flowRate_());
if (profile_.valid())
{
writeEntry(os, profile_());
}
if (!volumetric_)
{
writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
writeEntryIfDifferent<scalar>(os, "rhoInlet", -vGreat, rhoInlet_);
}
writeEntryIfDifferent<word>(os, "alpha", word::null, alphaName_);
writeEntry(os, "value", *this);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchVectorField,
flowRateInletVelocityFvPatchVectorField
);
}
// ************************************************************************* //