plenumPressureFvPatchScalarField: New plenum pressure boundary condition

This condition creates a zero-dimensional model of an enclosed volume of
gas upstream of the inlet. The pressure that the boundary condition
exerts on the inlet boundary is dependent on the thermodynamic state of
the upstream volume.  The upstream plenum density and temperature are
time-stepped along with the rest of the simulation, and momentum is
neglected. The plenum is supplied with a user specified mass flow and
temperature.

The result is a boundary condition which blends between a pressure inlet
condition condition and a fixed mass flow. The smaller the plenum
volume, the quicker the pressure responds to a deviation from the supply
mass flow, and the closer the model approximates a fixed mass flow. As
the plenum size increases, the model becomes more similar to a specified
pressure.

The expansion from the plenum to the inlet boundary is controlled by an
area ratio and a discharge coefficient. The area ratio can be used to
represent further acceleration between a sub-grid blockage such as fins.
The discharge coefficient represents a fractional deviation from an
ideal expansion process.

This condition is useful for simulating unsteady internal flow problems
for which both a mass flow boundary is unrealistic, and a pressure
boundary is susceptible to flow reversal. It was developed for use in
simulating confined combustion.

tutorials/compressible/rhoPimpleFoam/laminar/helmholtzResonance:
    helmholtz resonance tutorial case for plenum pressure boundary

This development was contributed by Will Bainbridge
This commit is contained in:
Henry Weller 2016-04-23 13:43:49 +01:00
parent 673e0d1704
commit 88561eea95
19 changed files with 1439 additions and 0 deletions

View File

@ -204,6 +204,7 @@ $(derivedFvPatchFields)/prghPressure/prghPressureFvPatchScalarField.C
$(derivedFvPatchFields)/prghTotalPressure/prghTotalPressureFvPatchScalarField.C
$(derivedFvPatchFields)/prghTotalHydrostaticPressure/prghTotalHydrostaticPressureFvPatchScalarField.C
$(derivedFvPatchFields)/fixedProfile/fixedProfileFvPatchFields.C
$(derivedFvPatchFields)/plenumPressure/plenumPressureFvPatchScalarField.C
fvsPatchFields = fields/fvsPatchFields
$(fvsPatchFields)/fvsPatchField/fvsPatchFields.C

View File

@ -0,0 +1,346 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 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 "plenumPressureFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
#include "surfaceFields.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::plenumPressureFvPatchScalarField::plenumPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchScalarField(p, iF),
gamma_(1.4),
R_(287.04),
supplyMassFlowRate_(1.0),
supplyTotalTemperature_(300.0),
plenumVolume_(1.0),
plenumDensity_(1.0),
plenumDensityOld_(1.0),
plenumTemperature_(300.0),
plenumTemperatureOld_(300.0),
rho_(1.0),
hasRho_(false),
inletAreaRatio_(1.0),
inletDischargeCoefficient_(1.0),
timeScale_(0.0),
timeIndex_(-1),
phiName_("phi"),
UName_("U")
{}
Foam::plenumPressureFvPatchScalarField::plenumPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const dictionary& dict
)
:
fixedValueFvPatchScalarField(p, iF),
gamma_(readScalar(dict.lookup("gamma"))),
R_(readScalar(dict.lookup("R"))),
supplyMassFlowRate_(readScalar(dict.lookup("supplyMassFlowRate"))),
supplyTotalTemperature_
(
readScalar(dict.lookup("supplyTotalTemperature"))
),
plenumVolume_(readScalar(dict.lookup("plenumVolume"))),
plenumDensity_(readScalar(dict.lookup("plenumDensity"))),
plenumTemperature_(readScalar(dict.lookup("plenumTemperature"))),
rho_(1.0),
hasRho_(false),
inletAreaRatio_(readScalar(dict.lookup("inletAreaRatio"))),
inletDischargeCoefficient_
(
readScalar(dict.lookup("inletDischargeCoefficient"))
),
timeScale_(dict.lookupOrDefault<scalar>("timeScale", 0.0)),
phiName_(dict.lookupOrDefault<word>("phi", "phi")),
UName_(dict.lookupOrDefault<word>("U", "U"))
{
fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
if (dict.found("rho"))
{
rho_ = readScalar(dict.lookup("rho"));
hasRho_ = true;
}
}
Foam::plenumPressureFvPatchScalarField::plenumPressureFvPatchScalarField
(
const plenumPressureFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
)
:
fixedValueFvPatchScalarField(ptf, p, iF, mapper),
gamma_(ptf.gamma_),
R_(ptf.R_),
supplyMassFlowRate_(ptf.supplyMassFlowRate_),
supplyTotalTemperature_(ptf.supplyTotalTemperature_),
plenumVolume_(ptf.plenumVolume_),
plenumDensity_(ptf.plenumDensity_),
plenumTemperature_(ptf.plenumTemperature_),
rho_(ptf.rho_),
hasRho_(ptf.hasRho_),
inletAreaRatio_(ptf.inletAreaRatio_),
inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
timeScale_(ptf.timeScale_),
phiName_(ptf.phiName_),
UName_(ptf.UName_)
{}
Foam::plenumPressureFvPatchScalarField::plenumPressureFvPatchScalarField
(
const plenumPressureFvPatchScalarField& tppsf
)
:
fixedValueFvPatchScalarField(tppsf),
gamma_(tppsf.gamma_),
R_(tppsf.R_),
supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
plenumVolume_(tppsf.plenumVolume_),
plenumDensity_(tppsf.plenumDensity_),
plenumTemperature_(tppsf.plenumTemperature_),
rho_(tppsf.rho_),
hasRho_(tppsf.hasRho_),
inletAreaRatio_(tppsf.inletAreaRatio_),
inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
timeScale_(tppsf.timeScale_),
phiName_(tppsf.phiName_),
UName_(tppsf.UName_)
{}
Foam::plenumPressureFvPatchScalarField::plenumPressureFvPatchScalarField
(
const plenumPressureFvPatchScalarField& tppsf,
const DimensionedField<scalar, volMesh>& iF
)
:
fixedValueFvPatchScalarField(tppsf, iF),
gamma_(tppsf.gamma_),
R_(tppsf.R_),
supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
plenumVolume_(tppsf.plenumVolume_),
plenumDensity_(tppsf.plenumDensity_),
plenumTemperature_(tppsf.plenumTemperature_),
rho_(tppsf.rho_),
hasRho_(tppsf.hasRho_),
inletAreaRatio_(tppsf.inletAreaRatio_),
inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
timeScale_(tppsf.timeScale_),
phiName_(tppsf.phiName_),
UName_(tppsf.UName_)
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::plenumPressureFvPatchScalarField::updateCoeffs()
{
if (updated())
{
return;
}
// Patch properties
const fvPatchField<scalar>& p = *this;
const fvPatchField<scalar>& p_old =
db().lookupObject<volScalarField>
(
dimensionedInternalField().name()
).oldTime().boundaryField()[patch().index()];
const fvPatchField<vector>& U =
patch().lookupPatchField<volVectorField, vector>(UName_);
const fvsPatchField<scalar>& phi =
patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
// Get the timestep
const scalar dt = db().time().deltaTValue();
// Check if operating at a new time index and update the old-time properties
// if so
if (timeIndex_ != db().time().timeIndex())
{
timeIndex_ = db().time().timeIndex();
plenumDensityOld_ = plenumDensity_;
plenumTemperatureOld_ = plenumTemperature_;
}
// Calculate the current mass flow rate
scalar massFlowRate(1.0);
if (phi.dimensionedInternalField().dimensions() == dimVelocity*dimArea)
{
if (hasRho_)
{
massFlowRate = - gSum(rho_*phi);
}
else
{
FatalErrorInFunction
<< "The density must be specified when using a volumetric flux."
<< exit(FatalError);
}
}
else if
(
phi.dimensionedInternalField().dimensions()
== dimDensity*dimVelocity*dimArea
)
{
if (hasRho_)
{
FatalErrorInFunction
<< "The density must be not specified when using a mass flux."
<< exit(FatalError);
}
else
{
massFlowRate = - gSum(phi);
}
}
else
{
FatalErrorInFunction
<< "dimensions of phi are not correct"
<< "\n on patch " << patch().name()
<< " of field " << dimensionedInternalField().name()
<< " in file " << dimensionedInternalField().objectPath() << nl
<< exit(FatalError);
}
// Calcaulate the specific heats
const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
// Calculate the new plenum properties
plenumDensity_ =
plenumDensityOld_
+ (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
plenumTemperature_ =
plenumTemperatureOld_
+ (dt/(plenumDensity_*cv*plenumVolume_))
* (
supplyMassFlowRate_
*(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
- massFlowRate*R_*plenumTemperature_
);
const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
// Squared velocity magnitude at exit of channels
const scalarField U_e(magSqr(U/inletAreaRatio_));
// Exit temperature to plenum temperature ratio
const scalarField r
(
1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
);
// Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
const scalarField a
(
(1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
);
// Isentropic exit temperature to plenum temperature ratio
const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
// Exit pressure to plenum pressure ratio
const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
// Limit to prevent outflow
const scalarField p_new
(
(1.0 - pos(phi))*t*plenumPressure + pos(phi)*max(p, plenumPressure)
);
// Relaxation fraction
const scalar oneByFraction = timeScale_/dt;
const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
// Set the new value
operator==((1.0 - fraction)*p_old + fraction*p_new);
fixedValueFvPatchScalarField::updateCoeffs();
}
void Foam::plenumPressureFvPatchScalarField::write(Ostream& os) const
{
fvPatchScalarField::write(os);
os.writeKeyword("gamma") << gamma_
<< token::END_STATEMENT << nl;
os.writeKeyword("R") << R_
<< token::END_STATEMENT << nl;
os.writeKeyword("supplyMassFlowRate") << supplyMassFlowRate_
<< token::END_STATEMENT << nl;
os.writeKeyword("supplyTotalTemperature") << supplyTotalTemperature_
<< token::END_STATEMENT << nl;
os.writeKeyword("plenumVolume") << plenumVolume_
<< token::END_STATEMENT << nl;
os.writeKeyword("plenumDensity") << plenumDensity_
<< token::END_STATEMENT << nl;
os.writeKeyword("plenumTemperature") << plenumTemperature_
<< token::END_STATEMENT << nl;
if (hasRho_)
{
os.writeKeyword("rho") << rho_
<< token::END_STATEMENT << nl;
}
os.writeKeyword("inletAreaRatio") << inletAreaRatio_
<< token::END_STATEMENT << nl;
os.writeKeyword("inletDischargeCoefficient") << inletDischargeCoefficient_
<< token::END_STATEMENT << nl;
writeEntryIfDifferent<scalar>(os, "timeScale", 0.0, timeScale_);
writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
writeEntryIfDifferent<word>(os, "U", "U", UName_);
writeEntry("value", os);
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
makePatchTypeField
(
fvPatchScalarField,
plenumPressureFvPatchScalarField
);
}
// ************************************************************************* //

View File

@ -0,0 +1,271 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016 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::plenumPressureFvPatchScalarField
Group
grpInletBoundaryConditions
Description
This boundary condition provides a plenum pressure inlet condition. This
condition creates a zero-dimensional model of an enclosed volume of gas
upstream of the inlet. The pressure that the boundary condition exerts on
the inlet boundary is dependent on the thermodynamic state of the upstream
volume. The upstream plenum density and temperature are time-stepped along
with the rest of the simulation, and momentum is neglected. The plenum is
supplied with a user specified mass flow and temperature.
The result is a boundary condition which blends between a pressure inlet
condition condition and a fixed mass flow. The smaller the plenum
volume, the quicker the pressure responds to a deviation from the supply
mass flow, and the closer the model approximates a fixed mass flow. As
the plenum size increases, the model becomes more similar to a specified
pressure.
The expansion from the plenum to the inlet boundary is controlled by an
area ratio and a discharge coefficient. The area ratio can be used to
represent further acceleration between a sub-grid blockage such as fins.
The discharge coefficient represents a fractional deviation from an
ideal expansion process.
This condition is useful for simulating unsteady internal flow problems
for which both a mass flow boundary is unrealistic, and a pressure
boundary is susceptible to flow reversal. It was developed for use in
simulating confined combustion.
Reference:
\verbatim
Bainbridge, W. (2013).
The Numerical Simulation of Oscillations in Gas Turbine Combustion
Chambers,
PhD Thesis,
Chapter 4, Section 4.3.1.2, 77-80.
\endverbatim
\heading Patch usage
\table
Property | Description | Required | Default value
gamma | ratio of specific heats | yes | none
R | specific gas constant | yes | none
supplyMassFlowRate | flow rate into the plenum | yes | none
supplyTotalTemperature | temperature into the plenum | yes | none
plenumVolume | plenum volume | yes | none
plenumDensity | plenum density | yes | none
plenumTemperature | plenum temperature | yes | none
U | velocity field name | no | U
phi | flux field name | no | phi
rho | inlet density | no | none
inletAreaRatio | inlet open fraction | yes | none
inletDischargeCoefficient | inlet loss coefficient | yes | none
timeScale | relaxation time scale | yes | none
\endtable
Example of the boundary condition specification:
\verbatim
myPatch
{
type plenumPressure;
gamma 1.4;
R 287.04;
supplyMassFlowRate 0.0001;
supplyTotalTemperature 300;
plenumVolume 0.000125;
plenumDensity 1.1613;
plenumTemperature 300;
inletAreaRatio 1.0;
inletDischargeCoefficient 0.8;
timeScale 1e-4;
value uniform 1e5;
}
\endverbatim
SourceFiles
plenumPressureFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef plenumPressureFvPatchScalarField_H
#define plenumPressureFvPatchScalarField_H
#include "fixedValueFvPatchFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class plenumPressureFvPatch Declaration
\*---------------------------------------------------------------------------*/
class plenumPressureFvPatchScalarField
:
public fixedValueFvPatchScalarField
{
// Private data
//- Ratio of specific heats
scalar gamma_;
//- Specific gas constant
scalar R_;
//- Mass flow rate supplied to the plenum
scalar supplyMassFlowRate_;
//- Total temperature of the gas supplied to the plenum
scalar supplyTotalTemperature_;
//- The volume of the plenum
scalar plenumVolume_;
//- The mean density of the gas in the plenum
scalar plenumDensity_;
//- The old-time mean density of the gas in the plenum
scalar plenumDensityOld_;
//- The mean temperature of the gas in the plenum
scalar plenumTemperature_;
//- The mean old-time temperature of the gas in the plenum
scalar plenumTemperatureOld_;
//- The constant density used when phi is volumetric
scalar rho_;
//- Whether or not the constant density has been specified
bool hasRho_;
//- The ratio of open area to total area at the inlet
// Allows a grid or mesh to be represented
scalar inletAreaRatio_;
//- The discharge coefficient at the inlet
scalar inletDischargeCoefficient_;
//- The time scale over which changes in pressure are smoothed
scalar timeScale_;
//- The time index used for updating
label timeIndex_;
//- The name of the flux field
word phiName_;
//- The name of the velocity field
word UName_;
public:
//- Runtime type information
TypeName("plenumPressure");
// Constructors
//- Construct from patch and internal field
plenumPressureFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
plenumPressureFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const dictionary&
);
//- Construct by mapping given plenumPressureFvPatchScalarField
// onto a new patch
plenumPressureFvPatchScalarField
(
const plenumPressureFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
plenumPressureFvPatchScalarField
(
const plenumPressureFvPatchScalarField&
);
//- Construct and return a clone
virtual tmp<fvPatchScalarField> clone() const
{
return tmp<fvPatchScalarField>
(
new plenumPressureFvPatchScalarField(*this)
);
}
//- Construct as copy setting internal field reference
plenumPressureFvPatchScalarField
(
const plenumPressureFvPatchScalarField&,
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 plenumPressureFvPatchScalarField(*this, iF)
);
}
// Member functions
// Evaluation functions
//- Update the coefficients associated with the patch field
virtual void updateCoeffs();
//- Write
virtual void write(Ostream&) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,55 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volScalarField;
location "0";
object T;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 0 0 1 0 0 0];
internalField uniform 300;
boundaryField
{
inlet
{
type fixedValue;
value uniform 300;
}
outlet
{
type inletOutlet;
inletValue uniform 300;
value uniform 300;
}
symmetry
{
type symmetry;
}
wall
{
type fixedValue;
value uniform 300;
}
plenum
{
type fixedValue;
value uniform 300;
}
}
// ************************************************************************* //

View File

@ -0,0 +1,60 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volVectorField;
location "0";
object U;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 1 -1 0 0 0 0];
internalField uniform (0 0 0);
boundaryField
{
inlet
{
type flowRateInletVelocity;
phi phi;
rho rho;
massFlowRate 0.0001;
value uniform (0 0 0);
}
outlet
{
type pressureInletOutletVelocity;
inletValue uniform (0 0 0);
value uniform (0 0 0);
}
symmetry
{
type symmetry;
}
wall
{
type fixedValue;
value uniform (0 0 0);
}
plenum
{
type pressureInletVelocity;
phi phi;
rho rho;
value uniform (0 0 0);
}
}
// ************************************************************************* //

View File

@ -0,0 +1,63 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volScalarField;
location "0";
object p;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [1 -1 -2 0 0 0 0];
internalField uniform 1e5;
boundaryField
{
inlet
{
type zeroGradient;
}
outlet
{
type fixedMean;
meanValue constant 1e5;
value uniform 1e5;
}
symmetry
{
type symmetry;
}
wall
{
type zeroGradient;
}
plenum
{
type plenumPressure;
gamma 1.4;
R 287.04;
supplyMassFlowRate 0.0001;
supplyTotalTemperature 300;
plenumVolume 0.000125;
plenumDensity 1.1613;
plenumTemperature 300;
inletAreaRatio 1.0;
inletDischargeCoefficient 0.8;
timeScale 1e-4;
value uniform 1e5;
}
}
// ************************************************************************* //

View File

@ -0,0 +1,9 @@
#!/bin/sh
. $WM_PROJECT_DIR/bin/tools/CleanFunctions
cleanCase
(cd system && rm -f blockMeshDict.caseBlocks blockMeshDict.caseBoundary)
rm -rf resolved modelled pressure.eps

View File

@ -0,0 +1,47 @@
#!/bin/sh
cd ${0%/*} || exit 1
. $WM_PROJECT_DIR/bin/tools/RunFunctions
# Run function links the appropriate mesh files and clones the case
run()
{
(
cd system
rm -f blockMeshDict.caseBlocks blockMeshDict.caseBoundary
ln -s blockMeshDict.${1}Blocks blockMeshDict.caseBlocks
ln -s blockMeshDict.${1}Boundary blockMeshDict.caseBoundary
)
cloneCase . ${1}
(
cd ${1}
runApplication blockMesh
runApplication decomposePar
runParallel $(getApplication)
)
}
# Run with a fully resolved plenum
run resolved
# Run with the plenum modelled by a boundary condition
run modelled
# Plot a comparison of the pressure in the neck
cat << EOF | gnuplot -persist
set terminal postscript eps size 5,4 enhanced color
set xlabel "Time (s)"
set ylabel "Guage pressure in the neck (Pa)"
set output "pressure.eps"
plot \
"resolved/postProcessing/probes/0/p" us 1:(\$4-1e5) t "Resolved Plenum" w l, \
"modelled/postProcessing/probes/0/p" us 1:(\$4-1e5) t "Modelled Plenum" w l
EOF

View File

@ -0,0 +1,49 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "constant";
object thermophysicalProperties;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
thermoType
{
type hePsiThermo;
mixture pureMixture;
transport const;
thermo eConst;
equationOfState perfectGas;
specie specie;
energy sensibleInternalEnergy;
}
mixture
{
specie
{
nMoles 1;
molWeight 28.9;
}
thermodynamics
{
Cv 712;
Hf 0;
}
transport
{
mu 1.8e-05;
Pr 0.7;
}
}
// ************************************************************************* //

View File

@ -0,0 +1,21 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "constant";
object turbulenceProperties;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
simulationType laminar;
// ************************************************************************* //

View File

@ -0,0 +1,251 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object blockMeshDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
convertToMeters 0.0025;
vertices
(
(-28 -10 -10)
(-28 -10 -5)
(-28 -10 5)
(-28 -10 10)
(-28 -5 -10)
(-28 -5 -5)
(-28 -5 5)
(-28 -5 10)
(-28 5 -10)
(-28 5 -5)
(-28 5 5)
(-28 5 10)
(-28 10 -10)
(-28 10 -5)
(-28 10 5)
(-28 10 10)
(-14 -10 -10)
(-14 -10 -5)
(-14 -10 5)
(-14 -10 10)
(-14 -5 -10)
(-14 -5 -5)
(-14 -5 5)
(-14 -5 10)
(-14 5 -10)
(-14 5 -5)
(-14 5 5)
(-14 5 10)
(-14 10 -10)
(-14 10 -5)
(-14 10 5)
(-14 10 10)
( -8 -10 -10)
( -8 -10 -5)
( -8 -10 5)
( -8 -10 10)
( -8 -5 -10)
( -8 -5 -5)
( -8 -5 5)
( -8 -5 10)
( -8 5 -10)
( -8 5 -5)
( -8 5 5)
( -8 5 10)
( -8 10 -10)
( -8 10 -5)
( -8 10 5)
( -8 10 10)
( 0 -10 -10)
( 0 -10 -5)
( 0 -10 5)
( 0 -10 10)
( 0 -5 -10)
( 0 -5 -5)
( 0 -5 5)
( 0 -5 10)
( 0 5 -10)
( 0 5 -5)
( 0 5 5)
( 0 5 10)
( 0 10 -10)
( 0 10 -5)
( 0 10 5)
( 0 10 10)
( 6 -10 -10)
( 6 -10 -5)
( 6 -10 5)
( 6 -10 10)
( 6 -5 -10)
( 6 -5 -5)
( 6 -5 5)
( 6 -5 10)
( 6 5 -10)
( 6 5 -5)
( 6 5 5)
( 6 5 10)
( 6 10 -10)
( 6 10 -5)
( 6 10 5)
( 6 10 10)
( 10 -10 -10)
( 10 -10 -5)
( 10 -10 5)
( 10 -10 10)
( 10 -5 -10)
( 10 -5 -5)
( 10 -5 5)
( 10 -5 10)
( 10 5 -10)
( 10 5 -5)
( 10 5 5)
( 10 5 10)
( 10 10 -10)
( 10 10 -5)
( 10 10 5)
( 10 10 10)
(-11 -2 -2)
(-11 -2 2)
(-11 2 -2)
(-11 2 2)
( -8 -2 -2)
( -8 -2 2)
( -8 2 -2)
( -8 2 2)
( 0 -2 -2)
( 0 -2 2)
( 0 2 -2)
( 0 2 2)
( 3 -2 -2)
( 3 -2 2)
( 3 2 -2)
( 3 2 2)
);
x1 14; x2 6; x3 20; x4 6; x5 4; // X divisions
yc 6; ys 12; zc $yc; zs $ys; // Y and Z corner and side divisions
o 6; // O-grid divisions
blocks
(
#include "blockMeshDict.caseBlocks"
hex ( 48 52 53 49 64 68 69 65) ($yc $zc $x4) simpleGrading (1 1 1)
hex ( 49 53 54 50 65 69 70 66) ($yc $zs $x4) simpleGrading (1 1 1)
hex ( 50 54 55 51 66 70 71 67) ($yc $zc $x4) simpleGrading (1 1 1)
hex ( 52 56 57 53 68 72 73 69) ($ys $zc $x4) simpleGrading (1 1 1)
hex ( 54 58 59 55 70 74 75 71) ($ys $zc $x4) simpleGrading (1 1 1)
hex ( 56 60 61 57 72 76 77 73) ($yc $zc $x4) simpleGrading (1 1 1)
hex ( 57 61 62 58 73 77 78 74) ($yc $zs $x4) simpleGrading (1 1 1)
hex ( 58 62 63 59 74 78 79 75) ($yc $zc $x4) simpleGrading (1 1 1)
hex ( 64 68 69 65 80 84 85 81) ($yc $zc $x5) simpleGrading (1 1 1)
hex ( 65 69 70 66 81 85 86 82) ($yc $zs $x5) simpleGrading (1 1 1)
hex ( 66 70 71 67 82 86 87 83) ($yc $zc $x5) simpleGrading (1 1 1)
hex ( 68 72 73 69 84 88 89 85) ($ys $zc $x5) simpleGrading (1 1 1)
hex ( 69 73 74 70 85 89 90 86) ($ys $zs $x5) simpleGrading (1 1 1)
hex ( 70 74 75 71 86 90 91 87) ($ys $zc $x5) simpleGrading (1 1 1)
hex ( 72 76 77 73 88 92 93 89) ($yc $zc $x5) simpleGrading (1 1 1)
hex ( 73 77 78 74 89 93 94 90) ($yc $zs $x5) simpleGrading (1 1 1)
hex ( 74 78 79 75 90 94 95 91) ($yc $zc $x5) simpleGrading (1 1 1)
hex (100 102 103 101 104 106 107 105) ($ys $zs $x3) simpleGrading (1 1 1)
hex (104 106 107 105 108 110 111 109) ($ys $zs $x4) simpleGrading (1 1 1)
hex (108 110 111 109 69 73 74 70 ) ($ys $zs $o ) simpleGrading (1 1 1)
hex ( 53 54 70 69 104 105 109 108) ($zs $x4 $o ) simpleGrading (1 1 1)
hex ( 57 53 69 73 106 104 108 110) ($ys $x4 $o ) simpleGrading (1 1 1)
hex ( 54 58 74 70 105 107 111 109) ($ys $x4 $o ) simpleGrading (1 1 1)
hex ( 58 57 73 74 107 106 110 111) ($zs $x4 $o ) simpleGrading (1 1 1)
);
edges
(
);
defaultPatch
{
name walls;
type wall;
}
boundary
(
#include "blockMeshDict.caseBoundary"
outlet
{
type patch;
faces
(
(80 84 85 81)
(81 85 86 82)
(82 86 87 83)
(84 88 89 85)
(85 89 90 86)
(86 90 91 87)
(88 92 93 89)
(89 93 94 90)
(90 94 95 91)
);
}
sides
{
type symmetry;
faces
(
(48 49 65 64)
(49 50 66 65)
(50 51 67 66)
(48 52 68 64)
(52 56 72 68)
(56 60 76 72)
(51 55 71 67)
(55 59 75 71)
(59 63 79 75)
(60 61 77 76)
(61 62 78 77)
(62 63 79 78)
(64 65 81 80)
(65 66 82 81)
(66 67 83 82)
(64 68 84 80)
(68 72 88 84)
(72 76 92 88)
(67 71 87 83)
(71 75 91 87)
(75 79 95 91)
(76 77 93 92)
(77 78 94 93)
(78 79 95 94)
);
}
);
mergePatchPairs
(
);
// ************************************************************************* //

View File

@ -0,0 +1,8 @@
plenum
{
type patch;
faces
(
(100 102 103 101)
);
}

View File

@ -0,0 +1,26 @@
hex ( 0 4 5 1 16 20 21 17) ($yc $zc $x1) simpleGrading (1 1 1)
hex ( 1 5 6 2 17 21 22 18) ($yc $zs $x1) simpleGrading (1 1 1)
hex ( 2 6 7 3 18 22 23 19) ($yc $zc $x1) simpleGrading (1 1 1)
hex ( 4 8 9 5 20 24 25 21) ($ys $zc $x1) simpleGrading (1 1 1)
hex ( 5 9 10 6 21 25 26 22) ($ys $zs $x1) simpleGrading (1 1 1)
hex ( 6 10 11 7 22 26 27 23) ($ys $zc $x1) simpleGrading (1 1 1)
hex ( 8 12 13 9 24 28 29 25) ($yc $zc $x1) simpleGrading (1 1 1)
hex ( 9 13 14 10 25 29 30 26) ($yc $zs $x1) simpleGrading (1 1 1)
hex ( 10 14 15 11 26 30 31 27) ($yc $zc $x1) simpleGrading (1 1 1)
hex ( 16 20 21 17 32 36 37 33) ($yc $zc $x2) simpleGrading (1 1 1)
hex ( 17 21 22 18 33 37 38 34) ($yc $zs $x2) simpleGrading (1 1 1)
hex ( 18 22 23 19 34 38 39 35) ($yc $zc $x2) simpleGrading (1 1 1)
hex ( 20 24 25 21 36 40 41 37) ($ys $zc $x2) simpleGrading (1 1 1)
hex ( 22 26 27 23 38 42 43 39) ($ys $zc $x2) simpleGrading (1 1 1)
hex ( 24 28 29 25 40 44 45 41) ($yc $zc $x2) simpleGrading (1 1 1)
hex ( 25 29 30 26 41 45 46 42) ($yc $zs $x2) simpleGrading (1 1 1)
hex ( 26 30 31 27 42 46 47 43) ($yc $zc $x2) simpleGrading (1 1 1)
hex ( 21 25 26 22 96 98 99 97) ($ys $zs $o ) simpleGrading (1 1 1)
hex ( 96 98 99 97 100 102 103 101) ($ys $zs $x2) simpleGrading (1 1 1)
hex ( 21 22 38 37 96 97 101 100) ($zs $x2 $o ) simpleGrading (1 1 1)
hex ( 25 21 37 41 98 96 100 102) ($ys $x2 $o ) simpleGrading (1 1 1)
hex ( 22 26 42 38 97 99 103 101) ($ys $x2 $o ) simpleGrading (1 1 1)
hex ( 26 25 41 42 99 98 102 103) ($zs $x2 $o ) simpleGrading (1 1 1)

View File

@ -0,0 +1,16 @@
inlet
{
type patch;
faces
(
( 0 4 5 1)
( 1 5 6 2)
( 2 6 7 3)
( 4 8 9 5)
( 5 9 10 6)
( 6 10 11 7)
( 8 12 13 9)
( 9 13 14 10)
(10 14 15 11)
);
}

View File

@ -0,0 +1,73 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object controlDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
application rhoPimpleFoam;
startFrom startTime;
startTime 0;
stopAt endTime;
endTime 0.05;
deltaT 1e-4;
writeControl adjustableRunTime;
writeInterval 1e-2;
purgeWrite 0;
writeFormat ascii;
writePrecision 10;
writeCompression off;
timeFormat general;
timePrecision 6;
runTimeModifiable true;
adjustTimeStep no;
maxCo 0.5;
functions
{
probes
{
functionObjectLibs ( "libsampling.so" );
type probes;
name probes;
outputControl timeStep;
outputInterval 1;
fields ( p );
probeLocations
(
( -0.045 0 0 )
( -0.045 0.020 0 )
( -0.010 0 0 )
( 0.0125 0 0 )
( 0.0125 0.020 0 )
);
}
}
// ************************************************************************* //

View File

@ -0,0 +1,22 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object decomposeParDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
numberOfSubdomains 4;
method scotch;
// ************************************************************************* //

View File

@ -0,0 +1,54 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object fvSchemes;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
ddtSchemes
{
default Euler;
}
gradSchemes
{
default Gauss linear;
}
divSchemes
{
default none;
div(phi,U) Gauss limitedLinearV 1;
div(phi,e) Gauss limitedLinear 1;
div(phi,K) Gauss limitedLinear 1;
div(phiv,p) Gauss limitedLinear 1;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
}
laplacianSchemes
{
default Gauss linear corrected;
}
interpolationSchemes
{
default linear;
}
snGradSchemes
{
default corrected;
}
// ************************************************************************* //

View File

@ -0,0 +1,67 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: dev |
| \\ / A nd | Web: www.OpenFOAM.org |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object fvSolution;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
solvers
{
"(p|rho)"
{
solver PCG;
preconditioner DIC;
tolerance 1e-6;
relTol 0.01;
}
"(p|rho)Final"
{
$p;
relTol 0;
}
"(U|e|k|nuTilda)"
{
solver smoothSolver;
smoother symGaussSeidel;
tolerance 1e-6;
relTol 0.01;
}
"(U|e|k|nuTilda)Final"
{
$U;
relTol 0;
}
}
PIMPLE
{
momentumPredictor yes;
nOuterCorrectors 3;
nCorrectors 1;
nNonOrthogonalCorrectors 0;
rhoMin 0.5;
rhoMax 2.0;
}
relaxationFactors
{
equations
{
".*" 1;
}
}
// ************************************************************************* //