ENH: Converted BC from timeVaryingUniformTotalPressure to uniformTotalPressure

This commit is contained in:
andy
2011-11-25 14:11:44 +00:00
parent ea52a8bca6
commit 169f15a593
3 changed files with 49 additions and 58 deletions

View File

@ -159,11 +159,11 @@ $(derivedFvPatchFields)/syringePressure/syringePressureFvPatchScalarField.C
$(derivedFvPatchFields)/timeVaryingMappedFixedValue/AverageIOFields.C
$(derivedFvPatchFields)/timeVaryingMappedFixedValue/timeVaryingMappedFixedValueFvPatchFields.C
$(derivedFvPatchFields)/totalPressure/totalPressureFvPatchScalarField.C
$(derivedFvPatchFields)/timeVaryingUniformTotalPressure/timeVaryingUniformTotalPressureFvPatchScalarField.C
$(derivedFvPatchFields)/totalTemperature/totalTemperatureFvPatchScalarField.C
$(derivedFvPatchFields)/turbulentInlet/turbulentInletFvPatchFields.C
$(derivedFvPatchFields)/turbulentIntensityKineticEnergyInlet/turbulentIntensityKineticEnergyInletFvPatchScalarField.C
$(derivedFvPatchFields)/uniformFixedValue/uniformFixedValueFvPatchFields.C
$(derivedFvPatchFields)/uniformTotalPressure/uniformTotalPressureFvPatchScalarField.C
$(derivedFvPatchFields)/waveTransmissive/waveTransmissiveFvPatchFields.C
$(derivedFvPatchFields)/uniformDensityHydrostaticPressure/uniformDensityHydrostaticPressureFvPatchScalarField.C
$(derivedFvPatchFields)/swirlFlowRateInletVelocity/swirlFlowRateInletVelocityFvPatchVectorField.C

View File

@ -23,7 +23,7 @@ License
\*---------------------------------------------------------------------------*/
#include "timeVaryingUniformTotalPressureFvPatchScalarField.H"
#include "uniformTotalPressureFvPatchScalarField.H"
#include "addToRunTimeSelectionTable.H"
#include "fvPatchFieldMapper.H"
#include "volFields.H"
@ -31,8 +31,8 @@ License
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::timeVaryingUniformTotalPressureFvPatchScalarField::
timeVaryingUniformTotalPressureFvPatchScalarField
Foam::uniformTotalPressureFvPatchScalarField::
uniformTotalPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF
@ -45,12 +45,12 @@ timeVaryingUniformTotalPressureFvPatchScalarField
psiName_("none"),
gamma_(0.0),
p0_(0.0),
timeSeries_()
pressure_()
{}
Foam::timeVaryingUniformTotalPressureFvPatchScalarField::
timeVaryingUniformTotalPressureFvPatchScalarField
Foam::uniformTotalPressureFvPatchScalarField::
uniformTotalPressureFvPatchScalarField
(
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
@ -64,7 +64,7 @@ timeVaryingUniformTotalPressureFvPatchScalarField
psiName_(dict.lookupOrDefault<word>("psi", "none")),
gamma_(readScalar(dict.lookup("gamma"))),
p0_(readScalar(dict.lookup("p0"))),
timeSeries_(dict)
pressure_(DataEntry<scalar>::New("pressure", dict))
{
if (dict.found("value"))
{
@ -80,10 +80,10 @@ timeVaryingUniformTotalPressureFvPatchScalarField
}
Foam::timeVaryingUniformTotalPressureFvPatchScalarField::
timeVaryingUniformTotalPressureFvPatchScalarField
Foam::uniformTotalPressureFvPatchScalarField::
uniformTotalPressureFvPatchScalarField
(
const timeVaryingUniformTotalPressureFvPatchScalarField& ptf,
const uniformTotalPressureFvPatchScalarField& ptf,
const fvPatch& p,
const DimensionedField<scalar, volMesh>& iF,
const fvPatchFieldMapper& mapper
@ -96,14 +96,14 @@ timeVaryingUniformTotalPressureFvPatchScalarField
psiName_(ptf.psiName_),
gamma_(ptf.gamma_),
p0_(ptf.p0_),
timeSeries_(ptf.timeSeries_)
pressure_(ptf.pressure_().clone().ptr())
{}
Foam::timeVaryingUniformTotalPressureFvPatchScalarField::
timeVaryingUniformTotalPressureFvPatchScalarField
Foam::uniformTotalPressureFvPatchScalarField::
uniformTotalPressureFvPatchScalarField
(
const timeVaryingUniformTotalPressureFvPatchScalarField& tppsf
const uniformTotalPressureFvPatchScalarField& tppsf
)
:
fixedValueFvPatchScalarField(tppsf),
@ -113,14 +113,14 @@ timeVaryingUniformTotalPressureFvPatchScalarField
psiName_(tppsf.psiName_),
gamma_(tppsf.gamma_),
p0_(tppsf.p0_),
timeSeries_(tppsf.timeSeries_)
pressure_(tppsf.pressure_().clone().ptr())
{}
Foam::timeVaryingUniformTotalPressureFvPatchScalarField::
timeVaryingUniformTotalPressureFvPatchScalarField
Foam::uniformTotalPressureFvPatchScalarField::
uniformTotalPressureFvPatchScalarField
(
const timeVaryingUniformTotalPressureFvPatchScalarField& tppsf,
const uniformTotalPressureFvPatchScalarField& tppsf,
const DimensionedField<scalar, volMesh>& iF
)
:
@ -131,13 +131,13 @@ timeVaryingUniformTotalPressureFvPatchScalarField
psiName_(tppsf.psiName_),
gamma_(tppsf.gamma_),
p0_(tppsf.p0_),
timeSeries_(tppsf.timeSeries_)
pressure_(tppsf.pressure_().clone().ptr())
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::timeVaryingUniformTotalPressureFvPatchScalarField::updateCoeffs
void Foam::uniformTotalPressureFvPatchScalarField::updateCoeffs
(
const vectorField& Up
)
@ -147,7 +147,7 @@ void Foam::timeVaryingUniformTotalPressureFvPatchScalarField::updateCoeffs
return;
}
p0_ = timeSeries_(this->db().time().timeOutputValue());
p0_ = pressure_->value(this->db().time().timeOutputValue());
const fvsPatchField<scalar>& phip =
patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
@ -189,10 +189,8 @@ void Foam::timeVaryingUniformTotalPressureFvPatchScalarField::updateCoeffs
}
else
{
FatalErrorIn
(
"timeVaryingUniformTotalPressureFvPatchScalarField::updateCoeffs()"
) << " rho or psi set inconsitently, rho = " << rhoName_
FatalErrorIn("uniformTotalPressureFvPatchScalarField::updateCoeffs()")
<< " rho or psi set inconsitently, rho = " << rhoName_
<< ", psi = " << psiName_ << ".\n"
<< " Set either rho or psi or neither depending on the "
"definition of total pressure.\n"
@ -207,14 +205,13 @@ void Foam::timeVaryingUniformTotalPressureFvPatchScalarField::updateCoeffs
}
void Foam::timeVaryingUniformTotalPressureFvPatchScalarField::updateCoeffs()
void Foam::uniformTotalPressureFvPatchScalarField::updateCoeffs()
{
updateCoeffs(patch().lookupPatchField<volVectorField, vector>(UName_));
}
void Foam::timeVaryingUniformTotalPressureFvPatchScalarField::
write(Ostream& os) const
void Foam::uniformTotalPressureFvPatchScalarField::write(Ostream& os) const
{
fvPatchScalarField::write(os);
writeEntryIfDifferent<word>(os, "U", "U", UName_);
@ -223,7 +220,7 @@ write(Ostream& os) const
os.writeKeyword("psi") << psiName_ << token::END_STATEMENT << nl;
os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
os.writeKeyword("p0") << p0_ << token::END_STATEMENT << nl;
timeSeries_.write(os);
pressure_->writeData(os);
writeEntry("value", os);
}
@ -235,7 +232,7 @@ namespace Foam
makePatchTypeField
(
fvPatchScalarField,
timeVaryingUniformTotalPressureFvPatchScalarField
uniformTotalPressureFvPatchScalarField
);
}

View File

@ -22,27 +22,27 @@ License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::timeVaryingUniformTotalPressureFvPatchScalarField
Foam::uniformTotalPressureFvPatchScalarField
Description
A time-varying form of a uniform total pressure boundary condition. The
variation is specified as an interpolationTable (see
Foam::interpolationTable).
variation is specified as an DataEntry (see Foam::DataEntry).
See Also
Foam::timeVaryingUniformFixedValueFvPatchField
Foam::uniformFixedValueFvPatchField.H
and Foam::totalPressureFvPatchScalarField.H
SourceFiles
timeVaryingUniformTotalPressureFvPatchScalarField.C
uniformTotalPressureFvPatchScalarField.C
\*---------------------------------------------------------------------------*/
#ifndef timeVaryingUniformTotalPressureFvPatchScalarField_H
#define timeVaryingUniformTotalPressureFvPatchScalarField_H
#ifndef uniformTotalPressureFvPatchScalarField_H
#define uniformTotalPressureFvPatchScalarField_H
#include "fixedValueFvPatchFields.H"
#include "interpolationTable.H"
#include "DataEntry.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -50,10 +50,10 @@ namespace Foam
{
/*---------------------------------------------------------------------------*\
Class timeVaryingTotalPressureFvPatch Declaration
Class uniformTotalPressureFvPatchField Declaration
\*---------------------------------------------------------------------------*/
class timeVaryingUniformTotalPressureFvPatchScalarField
class uniformTotalPressureFvPatchScalarField
:
public fixedValueFvPatchScalarField
{
@ -79,26 +79,26 @@ class timeVaryingUniformTotalPressureFvPatchScalarField
scalar p0_;
//- Table of time vs total pressure, including the bounding treatment
interpolationTable<scalar> timeSeries_;
autoPtr<DataEntry<scalar> > pressure_;
public:
//- Runtime type information
TypeName("timeVaryingTotalPressure");
TypeName("uniformTotalPressure");
// Constructors
//- Construct from patch and internal field
timeVaryingUniformTotalPressureFvPatchScalarField
uniformTotalPressureFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&
);
//- Construct from patch, internal field and dictionary
timeVaryingUniformTotalPressureFvPatchScalarField
uniformTotalPressureFvPatchScalarField
(
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
@ -106,18 +106,18 @@ public:
);
//- Construct by mapping given patch field onto a new patch
timeVaryingUniformTotalPressureFvPatchScalarField
uniformTotalPressureFvPatchScalarField
(
const timeVaryingUniformTotalPressureFvPatchScalarField&,
const uniformTotalPressureFvPatchScalarField&,
const fvPatch&,
const DimensionedField<scalar, volMesh>&,
const fvPatchFieldMapper&
);
//- Construct as copy
timeVaryingUniformTotalPressureFvPatchScalarField
uniformTotalPressureFvPatchScalarField
(
const timeVaryingUniformTotalPressureFvPatchScalarField&
const uniformTotalPressureFvPatchScalarField&
);
//- Construct and return a clone
@ -125,14 +125,14 @@ public:
{
return tmp<fvPatchScalarField>
(
new timeVaryingUniformTotalPressureFvPatchScalarField(*this)
new uniformTotalPressureFvPatchScalarField(*this)
);
}
//- Construct as copy setting internal field reference
timeVaryingUniformTotalPressureFvPatchScalarField
uniformTotalPressureFvPatchScalarField
(
const timeVaryingUniformTotalPressureFvPatchScalarField&,
const uniformTotalPressureFvPatchScalarField&,
const DimensionedField<scalar, volMesh>&
);
@ -144,7 +144,7 @@ public:
{
return tmp<fvPatchScalarField>
(
new timeVaryingUniformTotalPressureFvPatchScalarField(*this, iF)
new uniformTotalPressureFvPatchScalarField(*this, iF)
);
}
@ -190,12 +190,6 @@ public:
return p0_;
}
//- Return the time series used
const interpolationTable<scalar>& totalPressureTimeSeries() const
{
return timeSeries_;
}
// Evaluation functions