40 fixedValueFvPatchScalarField(p, iF),
43 supplyMassFlowRate_(1.0),
44 supplyTotalTemperature_(300.0),
47 plenumDensityOld_(1.0),
48 plenumTemperature_(300.0),
49 plenumTemperatureOld_(300.0),
53 inletDischargeCoefficient_(1.0),
68 fixedValueFvPatchScalarField(p, iF, dict),
69 gamma_(dict.
lookup<scalar>(
"gamma")),
70 R_(dict.
lookup<scalar>(
"R")),
71 supplyMassFlowRate_(dict.
lookup<scalar>(
"supplyMassFlowRate")),
72 supplyTotalTemperature_
74 dict.
lookup<scalar>(
"supplyTotalTemperature")
76 plenumVolume_(dict.
lookup<scalar>(
"plenumVolume")),
77 plenumDensity_(dict.
lookup<scalar>(
"plenumDensity")),
78 plenumTemperature_(dict.
lookup<scalar>(
"plenumTemperature")),
81 inletAreaRatio_(dict.
lookup<scalar>(
"inletAreaRatio")),
82 inletDischargeCoefficient_
84 dict.
lookup<scalar>(
"inletDischargeCoefficient")
90 if (dict.
found(
"rho"))
92 rho_ = dict.
lookup<scalar>(
"rho");
106 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
109 supplyMassFlowRate_(ptf.supplyMassFlowRate_),
110 supplyTotalTemperature_(ptf.supplyTotalTemperature_),
111 plenumVolume_(ptf.plenumVolume_),
112 plenumDensity_(ptf.plenumDensity_),
113 plenumTemperature_(ptf.plenumTemperature_),
115 hasRho_(ptf.hasRho_),
116 inletAreaRatio_(ptf.inletAreaRatio_),
117 inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
118 timeScale_(ptf.timeScale_),
119 phiName_(ptf.phiName_),
130 fixedValueFvPatchScalarField(tppsf, iF),
131 gamma_(tppsf.gamma_),
133 supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
134 supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
135 plenumVolume_(tppsf.plenumVolume_),
136 plenumDensity_(tppsf.plenumDensity_),
137 plenumTemperature_(tppsf.plenumTemperature_),
139 hasRho_(tppsf.hasRho_),
140 inletAreaRatio_(tppsf.inletAreaRatio_),
141 inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
142 timeScale_(tppsf.timeScale_),
143 phiName_(tppsf.phiName_),
162 internalField().name()
163 ).
oldTime().boundaryField()[patch().index()];
174 if (timeIndex_ != db().time().
timeIndex())
176 timeIndex_ = db().time().timeIndex();
177 plenumDensityOld_ = plenumDensity_;
178 plenumTemperatureOld_ = plenumTemperature_;
182 scalar massFlowRate(1.0);
183 if (phi.internalField().dimensions() ==
dimFlux)
187 massFlowRate = -
gSum(rho_*phi);
192 <<
"The density must be specified when using a volumetric flux." 198 phi.internalField().dimensions()
205 <<
"The density must be not specified when using a mass flux." 210 massFlowRate = -
gSum(phi);
216 <<
"dimensions of phi are not correct" 217 <<
"\n on patch " << patch().name()
218 <<
" of field " << internalField().name()
219 <<
" in file " << internalField().objectPath() <<
nl 224 const scalar cv = R_/(gamma_ - 1),
cp = R_*gamma_/(gamma_ - 1);
229 + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
231 plenumTemperatureOld_
232 + (dt/(plenumDensity_*cv*plenumVolume_))
235 *(
cp*supplyTotalTemperature_ - cv*plenumTemperature_)
236 - massFlowRate*R_*plenumTemperature_
238 const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
246 1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
252 (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
264 (1.0 -
pos0(phi))*t*plenumPressure +
pos0(phi)*
max(p, plenumPressure)
268 const scalar oneByFraction = timeScale_/dt;
269 const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
272 operator==((1.0 - fraction)*p_old + fraction*p_new);
273 fixedValueFvPatchScalarField::updateCoeffs();
282 writeEntry(os,
"supplyMassFlowRate", supplyMassFlowRate_);
283 writeEntry(os,
"supplyTotalTemperature", supplyTotalTemperature_);
284 writeEntry(os,
"plenumVolume", plenumVolume_);
285 writeEntry(os,
"plenumDensity", plenumDensity_);
286 writeEntry(os,
"plenumTemperature", plenumTemperature_);
291 writeEntry(os,
"inletAreaRatio", inletAreaRatio_);
295 "inletDischargeCoefficient",
296 inletDischargeCoefficient_
298 writeEntryIfDifferent<scalar>(os,
"timeScale", 0.0, timeScale_);
299 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
300 writeEntryIfDifferent<word>(os,
"U",
"U", UName_);
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
errorManipArg< error, int > exit(error &err, const int errNo=1)
A list of keyword definitions, which are a keyword followed by any number of values (e...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual void write(Ostream &) const
Write.
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy, recursively if necessary, the source to the destination.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Macros for easy insertion into run-time selection tables.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
scalar deltaTValue() const
Return time step value.
plenumPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
dimensionedScalar pos0(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const Time & time() const
Return time.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.