plenumPressureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2016-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  fixedValueFvPatchScalarField(p, iF),
41  gamma_(1.4),
42  R_(287.04),
43  supplyMassFlowRate_(1.0),
44  supplyTotalTemperature_(300.0),
45  plenumVolume_(1.0),
46  plenumDensity_(1.0),
47  plenumDensityOld_(1.0),
48  plenumTemperature_(300.0),
49  plenumTemperatureOld_(300.0),
50  rho_(1.0),
51  hasRho_(false),
52  inletAreaRatio_(1.0),
53  inletDischargeCoefficient_(1.0),
54  timeScale_(0.0),
55  timeIndex_(-1),
56  phiName_("phi"),
57  UName_("U")
58 {}
59 
60 
62 (
63  const fvPatch& p,
65  const dictionary& dict
66 )
67 :
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_
73  (
74  dict.lookup<scalar>("supplyTotalTemperature")
75  ),
76  plenumVolume_(dict.lookup<scalar>("plenumVolume")),
77  plenumDensity_(dict.lookup<scalar>("plenumDensity")),
78  plenumTemperature_(dict.lookup<scalar>("plenumTemperature")),
79  rho_(1.0),
80  hasRho_(false),
81  inletAreaRatio_(dict.lookup<scalar>("inletAreaRatio")),
82  inletDischargeCoefficient_
83  (
84  dict.lookup<scalar>("inletDischargeCoefficient")
85  ),
86  timeScale_(dict.lookupOrDefault<scalar>("timeScale", 0.0)),
87  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
88  UName_(dict.lookupOrDefault<word>("U", "U"))
89 {
90  if (dict.found("rho"))
91  {
92  rho_ = dict.lookup<scalar>("rho");
93  hasRho_ = true;
94  }
95 }
96 
97 
99 (
101  const fvPatch& p,
103  const fvPatchFieldMapper& mapper
104 )
105 :
106  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
107  gamma_(ptf.gamma_),
108  R_(ptf.R_),
109  supplyMassFlowRate_(ptf.supplyMassFlowRate_),
110  supplyTotalTemperature_(ptf.supplyTotalTemperature_),
111  plenumVolume_(ptf.plenumVolume_),
112  plenumDensity_(ptf.plenumDensity_),
113  plenumTemperature_(ptf.plenumTemperature_),
114  rho_(ptf.rho_),
115  hasRho_(ptf.hasRho_),
116  inletAreaRatio_(ptf.inletAreaRatio_),
117  inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
118  timeScale_(ptf.timeScale_),
119  phiName_(ptf.phiName_),
120  UName_(ptf.UName_)
121 {}
122 
123 
125 (
127 )
128 :
129  fixedValueFvPatchScalarField(tppsf),
130  gamma_(tppsf.gamma_),
131  R_(tppsf.R_),
132  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
133  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
134  plenumVolume_(tppsf.plenumVolume_),
135  plenumDensity_(tppsf.plenumDensity_),
136  plenumTemperature_(tppsf.plenumTemperature_),
137  rho_(tppsf.rho_),
138  hasRho_(tppsf.hasRho_),
139  inletAreaRatio_(tppsf.inletAreaRatio_),
140  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
141  timeScale_(tppsf.timeScale_),
142  phiName_(tppsf.phiName_),
143  UName_(tppsf.UName_)
144 {}
145 
146 
148 (
151 )
152 :
153  fixedValueFvPatchScalarField(tppsf, iF),
154  gamma_(tppsf.gamma_),
155  R_(tppsf.R_),
156  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
157  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
158  plenumVolume_(tppsf.plenumVolume_),
159  plenumDensity_(tppsf.plenumDensity_),
160  plenumTemperature_(tppsf.plenumTemperature_),
161  rho_(tppsf.rho_),
162  hasRho_(tppsf.hasRho_),
163  inletAreaRatio_(tppsf.inletAreaRatio_),
164  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
165  timeScale_(tppsf.timeScale_),
166  phiName_(tppsf.phiName_),
167  UName_(tppsf.UName_)
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
174 {
175  if (updated())
176  {
177  return;
178  }
179 
180  // Patch properties
181  const fvPatchField<scalar>& p = *this;
182  const fvPatchField<scalar>& p_old =
183  db().lookupObject<volScalarField>
184  (
185  internalField().name()
186  ).oldTime().boundaryField()[patch().index()];
187  const fvPatchField<vector>& U =
188  patch().lookupPatchField<volVectorField, vector>(UName_);
189  const fvsPatchField<scalar>& phi =
190  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
191 
192  // Get the timestep
193  const scalar dt = db().time().deltaTValue();
194 
195  // Check if operating at a new time index and update the old-time properties
196  // if so
197  if (timeIndex_ != db().time().timeIndex())
198  {
199  timeIndex_ = db().time().timeIndex();
200  plenumDensityOld_ = plenumDensity_;
201  plenumTemperatureOld_ = plenumTemperature_;
202  }
203 
204  // Calculate the current mass flow rate
205  scalar massFlowRate(1.0);
206  if (phi.internalField().dimensions() == dimVelocity*dimArea)
207  {
208  if (hasRho_)
209  {
210  massFlowRate = - gSum(rho_*phi);
211  }
212  else
213  {
215  << "The density must be specified when using a volumetric flux."
216  << exit(FatalError);
217  }
218  }
219  else if
220  (
221  phi.internalField().dimensions()
223  )
224  {
225  if (hasRho_)
226  {
228  << "The density must be not specified when using a mass flux."
229  << exit(FatalError);
230  }
231  else
232  {
233  massFlowRate = - gSum(phi);
234  }
235  }
236  else
237  {
239  << "dimensions of phi are not correct"
240  << "\n on patch " << patch().name()
241  << " of field " << internalField().name()
242  << " in file " << internalField().objectPath() << nl
243  << exit(FatalError);
244  }
245 
246  // Calcaulate the specific heats
247  const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
248 
249  // Calculate the new plenum properties
250  plenumDensity_ =
251  plenumDensityOld_
252  + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
253  plenumTemperature_ =
254  plenumTemperatureOld_
255  + (dt/(plenumDensity_*cv*plenumVolume_))
256  * (
257  supplyMassFlowRate_
258  *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
259  - massFlowRate*R_*plenumTemperature_
260  );
261  const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
262 
263  // Squared velocity magnitude at exit of channels
264  const scalarField U_e(magSqr(U/inletAreaRatio_));
265 
266  // Exit temperature to plenum temperature ratio
267  const scalarField r
268  (
269  1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
270  );
271 
272  // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
273  const scalarField a
274  (
275  (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
276  );
277 
278  // Isentropic exit temperature to plenum temperature ratio
279  const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
280 
281  // Exit pressure to plenum pressure ratio
282  const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
283 
284  // Limit to prevent outflow
285  const scalarField p_new
286  (
287  (1.0 - pos0(phi))*t*plenumPressure + pos0(phi)*max(p, plenumPressure)
288  );
289 
290  // Relaxation fraction
291  const scalar oneByFraction = timeScale_/dt;
292  const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
293 
294  // Set the new value
295  operator==((1.0 - fraction)*p_old + fraction*p_new);
296  fixedValueFvPatchScalarField::updateCoeffs();
297 }
298 
299 
301 {
303  writeEntry(os, "gamma", gamma_);
304  writeEntry(os, "R", R_);
305  writeEntry(os, "supplyMassFlowRate", supplyMassFlowRate_);
306  writeEntry(os, "supplyTotalTemperature", supplyTotalTemperature_);
307  writeEntry(os, "plenumVolume", plenumVolume_);
308  writeEntry(os, "plenumDensity", plenumDensity_);
309  writeEntry(os, "plenumTemperature", plenumTemperature_);
310  if (hasRho_)
311  {
312  writeEntry(os, "rho", rho_);
313  }
314  writeEntry(os, "inletAreaRatio", inletAreaRatio_);
315  writeEntry
316  (
317  os,
318  "inletDischargeCoefficient",
319  inletDischargeCoefficient_
320  );
321  writeEntryIfDifferent<scalar>(os, "timeScale", 0.0, timeScale_);
322  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
323  writeEntryIfDifferent<word>(os, "U", "U", UName_);
324  writeEntry(os, "value", *this);
325 }
326 
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 namespace Foam
331 {
333  (
336  );
337 }
338 
339 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
rho oldTime()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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.
Definition: POSIX.C:753
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Macros for easy insertion into run-time selection tables.
phi
Definition: pEqn.H:104
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.
Definition: word.H:59
Foam::fvPatchFieldMapper.
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
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...
Definition: Ostream.H:54
dimensionedScalar pos0(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:260
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const Time & time() const
Return time.
Definition: IOobject.C:328
label timeIndex
Definition: getTimeIndex.H:4
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...
Definition: fvsPatchField.H:65
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
const dimensionSet dimVelocity