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-2021 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 (
128 )
129 :
130  fixedValueFvPatchScalarField(tppsf, iF),
131  gamma_(tppsf.gamma_),
132  R_(tppsf.R_),
133  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
134  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
135  plenumVolume_(tppsf.plenumVolume_),
136  plenumDensity_(tppsf.plenumDensity_),
137  plenumTemperature_(tppsf.plenumTemperature_),
138  rho_(tppsf.rho_),
139  hasRho_(tppsf.hasRho_),
140  inletAreaRatio_(tppsf.inletAreaRatio_),
141  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
142  timeScale_(tppsf.timeScale_),
143  phiName_(tppsf.phiName_),
144  UName_(tppsf.UName_)
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 {
152  if (updated())
153  {
154  return;
155  }
156 
157  // Patch properties
158  const fvPatchField<scalar>& p = *this;
159  const fvPatchField<scalar>& p_old =
160  db().lookupObject<volScalarField>
161  (
162  internalField().name()
163  ).oldTime().boundaryField()[patch().index()];
164  const fvPatchField<vector>& U =
165  patch().lookupPatchField<volVectorField, vector>(UName_);
166  const fvsPatchField<scalar>& phi =
167  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
168 
169  // Get the timestep
170  const scalar dt = db().time().deltaTValue();
171 
172  // Check if operating at a new time index and update the old-time properties
173  // if so
174  if (timeIndex_ != db().time().timeIndex())
175  {
176  timeIndex_ = db().time().timeIndex();
177  plenumDensityOld_ = plenumDensity_;
178  plenumTemperatureOld_ = plenumTemperature_;
179  }
180 
181  // Calculate the current mass flow rate
182  scalar massFlowRate(1.0);
183  if (phi.internalField().dimensions() == dimFlux)
184  {
185  if (hasRho_)
186  {
187  massFlowRate = - gSum(rho_*phi);
188  }
189  else
190  {
192  << "The density must be specified when using a volumetric flux."
193  << exit(FatalError);
194  }
195  }
196  else if
197  (
198  phi.internalField().dimensions()
199  == dimMassFlux
200  )
201  {
202  if (hasRho_)
203  {
205  << "The density must be not specified when using a mass flux."
206  << exit(FatalError);
207  }
208  else
209  {
210  massFlowRate = - gSum(phi);
211  }
212  }
213  else
214  {
216  << "dimensions of phi are not correct"
217  << "\n on patch " << patch().name()
218  << " of field " << internalField().name()
219  << " in file " << internalField().objectPath() << nl
220  << exit(FatalError);
221  }
222 
223  // Calcaulate the specific heats
224  const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
225 
226  // Calculate the new plenum properties
227  plenumDensity_ =
228  plenumDensityOld_
229  + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
230  plenumTemperature_ =
231  plenumTemperatureOld_
232  + (dt/(plenumDensity_*cv*plenumVolume_))
233  * (
234  supplyMassFlowRate_
235  *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
236  - massFlowRate*R_*plenumTemperature_
237  );
238  const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
239 
240  // Squared velocity magnitude at exit of channels
241  const scalarField U_e(magSqr(U/inletAreaRatio_));
242 
243  // Exit temperature to plenum temperature ratio
244  const scalarField r
245  (
246  1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
247  );
248 
249  // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
250  const scalarField a
251  (
252  (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
253  );
254 
255  // Isentropic exit temperature to plenum temperature ratio
256  const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
257 
258  // Exit pressure to plenum pressure ratio
259  const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
260 
261  // Limit to prevent outflow
262  const scalarField p_new
263  (
264  (1.0 - pos0(phi))*t*plenumPressure + pos0(phi)*max(p, plenumPressure)
265  );
266 
267  // Relaxation fraction
268  const scalar oneByFraction = timeScale_/dt;
269  const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
270 
271  // Set the new value
272  operator==((1.0 - fraction)*p_old + fraction*p_new);
273  fixedValueFvPatchScalarField::updateCoeffs();
274 }
275 
276 
278 {
280  writeEntry(os, "gamma", gamma_);
281  writeEntry(os, "R", R_);
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_);
287  if (hasRho_)
288  {
289  writeEntry(os, "rho", rho_);
290  }
291  writeEntry(os, "inletAreaRatio", inletAreaRatio_);
292  writeEntry
293  (
294  os,
295  "inletDischargeCoefficient",
296  inletDischargeCoefficient_
297  );
298  writeEntryIfDifferent<scalar>(os, "timeScale", 0.0, timeScale_);
299  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
300  writeEntryIfDifferent<word>(os, "U", "U", UName_);
301  writeEntry(os, "value", *this);
302 }
303 
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 namespace Foam
308 {
310  (
313  );
314 }
315 
316 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
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)
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:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
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:62
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:260
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
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.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
phi
Definition: correctPhi.H:3
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
rho oldTime()
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimMassFlux
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:318
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
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844