Skip to content
Snippets Groups Projects
Commit e757cbd0 authored by sergio's avatar sergio
Browse files

ENH: swirlMassFlowRateInletVelocity inlet BC

parent bee7ba14
Branches
Tags
No related merge requests found
...@@ -150,6 +150,7 @@ $(derivedFvPatchFields)/turbulentIntensityKineticEnergyInlet/turbulentIntensityK ...@@ -150,6 +150,7 @@ $(derivedFvPatchFields)/turbulentIntensityKineticEnergyInlet/turbulentIntensityK
$(derivedFvPatchFields)/uniformFixedValue/uniformFixedValueFvPatchFields.C $(derivedFvPatchFields)/uniformFixedValue/uniformFixedValueFvPatchFields.C
$(derivedFvPatchFields)/waveTransmissive/waveTransmissiveFvPatchFields.C $(derivedFvPatchFields)/waveTransmissive/waveTransmissiveFvPatchFields.C
$(derivedFvPatchFields)/uniformDensityHydrostaticPressure/uniformDensityHydrostaticPressureFvPatchScalarField.C $(derivedFvPatchFields)/uniformDensityHydrostaticPressure/uniformDensityHydrostaticPressureFvPatchScalarField.C
$(derivedFvPatchFields)/swirlMassFlowRateInletVelocity/swirlMassFlowRateInletVelocityFvPatchVectorField.C
fvsPatchFields = fields/fvsPatchFields fvsPatchFields = fields/fvsPatchFields
$(fvsPatchFields)/fvsPatchField/fvsPatchFields.C $(fvsPatchFields)/fvsPatchField/fvsPatchFields.C
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2006-2009 OpenCFD Ltd.
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
\*---------------------------------------------------------------------------*/
#include "swirlMassFlowRateInletVelocityFvPatchVectorField.H"
#include "volFields.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "surfaceFields.H"
#include "mathematicalConstants.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::
swirlMassFlowRateInletVelocityFvPatchVectorField::
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF
)
:
fixedValueFvPatchField<vector>(p, iF),
flowRate_(0),
phiName_("phi"),
rhoName_("rho"),
rpm_(0)
{}
Foam::
swirlMassFlowRateInletVelocityFvPatchVectorField::
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const swirlMassFlowRateInletVelocityFvPatchVectorField& ptf,
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
flowRate_(ptf.flowRate_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_)
{}
Foam::
swirlMassFlowRateInletVelocityFvPatchVectorField::
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const fvPatch& p,
const DimensionedField<vector, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchField<vector>(p, iF, dict),
flowRate_(readScalar(dict.lookup("flowRate"))),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
rpm_(readScalar(dict.lookup("rpm")))
{}
Foam::
swirlMassFlowRateInletVelocityFvPatchVectorField::
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const swirlMassFlowRateInletVelocityFvPatchVectorField& ptf
)
:
fixedValueFvPatchField<vector>(ptf),
flowRate_(ptf.flowRate_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
rpm_(ptf.rpm_)
{}
Foam::
swirlMassFlowRateInletVelocityFvPatchVectorField::
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const swirlMassFlowRateInletVelocityFvPatchVectorField& ptf,
const DimensionedField<vector, volMesh>& iF
)
:
fixedValueFvPatchField<vector>(ptf, iF),
flowRate_(ptf.flowRate_),
phiName_(ptf.phiName_),
rhoName_(ptf.rhoName_),
rpm_(ptf.rpm_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::swirlMassFlowRateInletVelocityFvPatchVectorField::updateCoeffs()
{
if (updated())
{
return;
}
scalar totArea = gSum(patch().magSf());
// a simpler way of doing this would be nice
scalar avgU = -flowRate_/totArea;
vector center = gSum(patch().Cf()*patch().magSf())/totArea;
vector normal = gSum(patch().nf()*patch().magSf())/totArea;
vectorField tangVelo =
(rpm_*constant::mathematical::pi/30.0)
*(patch().Cf() - center) ^ normal;
vectorField n = patch().nf();
const surfaceScalarField& phi =
db().lookupObject<surfaceScalarField>(phiName_);
if (phi.dimensions() == dimVelocity*dimArea)
{
// volumetric flow-rate
operator==(tangVelo + n*avgU);
}
else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
{
const fvPatchField<scalar>& rhop =
patch().lookupPatchField<volScalarField, scalar>(rhoName_);
// mass flow-rate
operator==(tangVelo + n*avgU/rhop);
}
else
{
FatalErrorIn
(
"swirlMassFlowRateInletVelocityFvPatchVectorField::updateCoeffs()"
) << "dimensions of " << phiName_ << " are incorrect" << nl
<< " on patch " << this->patch().name()
<< " of field " << this->dimensionedInternalField().name()
<< " in file " << this->dimensionedInternalField().objectPath()
<< nl << exit(FatalError);
}
fixedValueFvPatchField<vector>::updateCoeffs();
}
void Foam::swirlMassFlowRateInletVelocityFvPatchVectorField::write(Ostream& os) const
{
fvPatchField<vector>::write(os);
os.writeKeyword("flowRate") << flowRate_ << token::END_STATEMENT << nl;
if (phiName_ != "phi")
{
os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
}
if (rhoName_ != "rho")
{
os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
}
os.writeKeyword("rpm") << rpm_ << token::END_STATEMENT << nl;
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchVectorField,
swirlMassFlowRateInletVelocityFvPatchVectorField
);
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2010-2010 OpenCFD Ltd.
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Class
Foam::swirlMassFlowRateInletVelocityFvPatchVectorField
Description
Describes a volumetric/mass flow normal vector boundary condition by its
magnitude as an integral over its area with a swirl component determined
by the RPM
The basis of the patch (volumetric or mass) is determined by the
dimensions of the flux, phi.
The current density is used to correct the velocity when applying the
mass basis.
Example of the boundary condition specification:
@verbatim
inlet
{
type swirlMassFlowRateInletVelocity;
flowRate 0.2; // Volumetric/mass flow rate [m3/s or kg/s]
rpm 100;
}
@endverbatim
Note
- The value is positive inwards
SourceFiles
swirlMassFlowRateInletVelocityFvPatchVectorField.C
\*---------------------------------------------------------------------------*/
#ifndef swirlMassFlowRateInletVelocityFvPatchVectorField_H
#define swirlMassFlowRateInletVelocityFvPatchVectorField_H
#include "fixedValueFvPatchFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class swirlMassFlowRateInletVelocityFvPatchVectorField Declaration
\*---------------------------------------------------------------------------*/
class swirlMassFlowRateInletVelocityFvPatchVectorField
:
public fixedValueFvPatchVectorField
{
// Private data
//- Inlet integral flow rate
scalar flowRate_;
//- Name of the flux transporting the field
word phiName_;
//- Name of the density field used to normalize the mass flux
word rhoName_;
//- RPM
scalar rpm_;
public:
//- Runtime type information
TypeName("swirlMassFlowRateInletVelocity");
// Constructors
//- Construct from patch and internal field
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const fvPatch&,
const DimensionedField<vector, volMesh>&
);
//- Construct from patch, internal field and dictionary
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const fvPatch&,
const DimensionedField<vector, volMesh>&,
const dictionary&
);
//- Construct by mapping given
// flowRateInletVelocityFvPatchVectorField
// onto a new patch
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const swirlMassFlowRateInletVelocityFvPatchVectorField&,
const fvPatch&,
const DimensionedField<vector, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const swirlMassFlowRateInletVelocityFvPatchVectorField&
);
//- Construct and return a clone
virtual tmp<fvPatchVectorField> clone() const
{
return tmp<fvPatchVectorField>
(
new swirlMassFlowRateInletVelocityFvPatchVectorField(*this)
);
}
//- Construct as copy setting internal field reference
swirlMassFlowRateInletVelocityFvPatchVectorField
(
const swirlMassFlowRateInletVelocityFvPatchVectorField&,
const DimensionedField<vector, volMesh>&
);
//- Construct and return a clone setting internal field reference
virtual tmp<fvPatchVectorField> clone
(
const DimensionedField<vector, volMesh>& iF
) const
{
return tmp<fvPatchVectorField>
(
new swirlMassFlowRateInletVelocityFvPatchVectorField(*this, iF)
);
}
// Member functions
// Access
//- Return the flux
scalar flowRate() const
{
return flowRate_;
}
//- Return reference to the flux to allow adjustment
scalar& flowRate()
{
return flowRate_;
}
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment