plenumPressureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2016 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),
69  gamma_(readScalar(dict.lookup("gamma"))),
70  R_(readScalar(dict.lookup("R"))),
71  supplyMassFlowRate_(readScalar(dict.lookup("supplyMassFlowRate"))),
72  supplyTotalTemperature_
73  (
74  readScalar(dict.lookup("supplyTotalTemperature"))
75  ),
76  plenumVolume_(readScalar(dict.lookup("plenumVolume"))),
77  plenumDensity_(readScalar(dict.lookup("plenumDensity"))),
78  plenumTemperature_(readScalar(dict.lookup("plenumTemperature"))),
79  rho_(1.0),
80  hasRho_(false),
81  inletAreaRatio_(readScalar(dict.lookup("inletAreaRatio"))),
82  inletDischargeCoefficient_
83  (
84  readScalar(dict.lookup("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  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
91 
92  if (dict.found("rho"))
93  {
94  rho_ = readScalar(dict.lookup("rho"));
95  hasRho_ = true;
96  }
97 }
98 
99 
101 (
103  const fvPatch& p,
105  const fvPatchFieldMapper& mapper
106 )
107 :
108  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
109  gamma_(ptf.gamma_),
110  R_(ptf.R_),
111  supplyMassFlowRate_(ptf.supplyMassFlowRate_),
112  supplyTotalTemperature_(ptf.supplyTotalTemperature_),
113  plenumVolume_(ptf.plenumVolume_),
114  plenumDensity_(ptf.plenumDensity_),
115  plenumTemperature_(ptf.plenumTemperature_),
116  rho_(ptf.rho_),
117  hasRho_(ptf.hasRho_),
118  inletAreaRatio_(ptf.inletAreaRatio_),
119  inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
120  timeScale_(ptf.timeScale_),
121  phiName_(ptf.phiName_),
122  UName_(ptf.UName_)
123 {}
124 
125 
127 (
129 )
130 :
131  fixedValueFvPatchScalarField(tppsf),
132  gamma_(tppsf.gamma_),
133  R_(tppsf.R_),
134  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
135  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
136  plenumVolume_(tppsf.plenumVolume_),
137  plenumDensity_(tppsf.plenumDensity_),
138  plenumTemperature_(tppsf.plenumTemperature_),
139  rho_(tppsf.rho_),
140  hasRho_(tppsf.hasRho_),
141  inletAreaRatio_(tppsf.inletAreaRatio_),
142  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
143  timeScale_(tppsf.timeScale_),
144  phiName_(tppsf.phiName_),
145  UName_(tppsf.UName_)
146 {}
147 
148 
150 (
153 )
154 :
155  fixedValueFvPatchScalarField(tppsf, iF),
156  gamma_(tppsf.gamma_),
157  R_(tppsf.R_),
158  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
159  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
160  plenumVolume_(tppsf.plenumVolume_),
161  plenumDensity_(tppsf.plenumDensity_),
162  plenumTemperature_(tppsf.plenumTemperature_),
163  rho_(tppsf.rho_),
164  hasRho_(tppsf.hasRho_),
165  inletAreaRatio_(tppsf.inletAreaRatio_),
166  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
167  timeScale_(tppsf.timeScale_),
168  phiName_(tppsf.phiName_),
169  UName_(tppsf.UName_)
170 {}
171 
172 
173 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
174 
176 {
177  if (updated())
178  {
179  return;
180  }
181 
182  // Patch properties
183  const fvPatchField<scalar>& p = *this;
184  const fvPatchField<scalar>& p_old =
185  db().lookupObject<volScalarField>
186  (
187  internalField().name()
188  ).oldTime().boundaryField()[patch().index()];
189  const fvPatchField<vector>& U =
190  patch().lookupPatchField<volVectorField, vector>(UName_);
191  const fvsPatchField<scalar>& phi =
192  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
193 
194  // Get the timestep
195  const scalar dt = db().time().deltaTValue();
196 
197  // Check if operating at a new time index and update the old-time properties
198  // if so
199  if (timeIndex_ != db().time().timeIndex())
200  {
201  timeIndex_ = db().time().timeIndex();
202  plenumDensityOld_ = plenumDensity_;
203  plenumTemperatureOld_ = plenumTemperature_;
204  }
205 
206  // Calculate the current mass flow rate
207  scalar massFlowRate(1.0);
208  if (phi.internalField().dimensions() == dimVelocity*dimArea)
209  {
210  if (hasRho_)
211  {
212  massFlowRate = - gSum(rho_*phi);
213  }
214  else
215  {
217  << "The density must be specified when using a volumetric flux."
218  << exit(FatalError);
219  }
220  }
221  else if
222  (
223  phi.internalField().dimensions()
225  )
226  {
227  if (hasRho_)
228  {
230  << "The density must be not specified when using a mass flux."
231  << exit(FatalError);
232  }
233  else
234  {
235  massFlowRate = - gSum(phi);
236  }
237  }
238  else
239  {
241  << "dimensions of phi are not correct"
242  << "\n on patch " << patch().name()
243  << " of field " << internalField().name()
244  << " in file " << internalField().objectPath() << nl
245  << exit(FatalError);
246  }
247 
248  // Calcaulate the specific heats
249  const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
250 
251  // Calculate the new plenum properties
252  plenumDensity_ =
253  plenumDensityOld_
254  + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
255  plenumTemperature_ =
256  plenumTemperatureOld_
257  + (dt/(plenumDensity_*cv*plenumVolume_))
258  * (
259  supplyMassFlowRate_
260  *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
261  - massFlowRate*R_*plenumTemperature_
262  );
263  const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
264 
265  // Squared velocity magnitude at exit of channels
266  const scalarField U_e(magSqr(U/inletAreaRatio_));
267 
268  // Exit temperature to plenum temperature ratio
269  const scalarField r
270  (
271  1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
272  );
273 
274  // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
275  const scalarField a
276  (
277  (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
278  );
279 
280  // Isentropic exit temperature to plenum temperature ratio
281  const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
282 
283  // Exit pressure to plenum pressure ratio
284  const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
285 
286  // Limit to prevent outflow
287  const scalarField p_new
288  (
289  (1.0 - pos(phi))*t*plenumPressure + pos(phi)*max(p, plenumPressure)
290  );
291 
292  // Relaxation fraction
293  const scalar oneByFraction = timeScale_/dt;
294  const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
295 
296  // Set the new value
297  operator==((1.0 - fraction)*p_old + fraction*p_new);
298  fixedValueFvPatchScalarField::updateCoeffs();
299 }
300 
301 
303 {
305  os.writeKeyword("gamma") << gamma_
306  << token::END_STATEMENT << nl;
307  os.writeKeyword("R") << R_
308  << token::END_STATEMENT << nl;
309  os.writeKeyword("supplyMassFlowRate") << supplyMassFlowRate_
310  << token::END_STATEMENT << nl;
311  os.writeKeyword("supplyTotalTemperature") << supplyTotalTemperature_
312  << token::END_STATEMENT << nl;
313  os.writeKeyword("plenumVolume") << plenumVolume_
314  << token::END_STATEMENT << nl;
315  os.writeKeyword("plenumDensity") << plenumDensity_
316  << token::END_STATEMENT << nl;
317  os.writeKeyword("plenumTemperature") << plenumTemperature_
318  << token::END_STATEMENT << nl;
319  if (hasRho_)
320  {
321  os.writeKeyword("rho") << rho_
322  << token::END_STATEMENT << nl;
323  }
324  os.writeKeyword("inletAreaRatio") << inletAreaRatio_
325  << token::END_STATEMENT << nl;
326  os.writeKeyword("inletDischargeCoefficient") << inletDischargeCoefficient_
327  << token::END_STATEMENT << nl;
328  writeEntryIfDifferent<scalar>(os, "timeScale", 0.0, timeScale_);
329  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
330  writeEntryIfDifferent<word>(os, "U", "U", UName_);
331  writeEntry("value", os);
332 }
333 
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 namespace Foam
338 {
340  (
343  );
344 }
345 
346 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
Foam::surfaceFields.
surfaceScalarField & phi
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
rho oldTime()
U
Definition: pEqn.H:83
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:137
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)
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:65
bool cp(const fileName &src, const fileName &dst)
Copy, recursively if necessary, the source to the destination.
Definition: POSIX.C:620
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dimensionedScalar pos(const dimensionedScalar &ds)
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
const Time & time() const
Return time.
Definition: IOobject.C:227
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
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:53
volScalarField scalarField(fieldObject, mesh)
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
const dimensionSet dimDensity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
const dimensionSet dimVelocity