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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF),
43  gamma_(1.4),
44  R_(287.04),
45  supplyMassFlowRate_(1.0),
46  supplyTotalTemperature_(300.0),
47  plenumVolume_(1.0),
48  plenumDensity_(1.0),
49  plenumDensityOld_(1.0),
50  plenumTemperature_(300.0),
51  plenumTemperatureOld_(300.0),
52  rho_(1.0),
53  hasRho_(false),
54  inletAreaRatio_(1.0),
55  inletDischargeCoefficient_(1.0),
56  timeScale_(0.0),
57  timeIndex_(-1),
58  phiName_("phi"),
59  UName_("U")
60 {}
61 
62 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  fixedValueFvPatchScalarField(p, iF, dict),
71  gamma_(dict.get<scalar>("gamma")),
72  R_(dict.get<scalar>("R")),
73  supplyMassFlowRate_(dict.get<scalar>("supplyMassFlowRate")),
74  supplyTotalTemperature_(dict.get<scalar>("supplyTotalTemperature")),
75  plenumVolume_(dict.get<scalar>("plenumVolume")),
76  plenumDensity_(dict.get<scalar>("plenumDensity")),
77  plenumTemperature_(dict.get<scalar>("plenumTemperature")),
78  rho_(1.0),
79  hasRho_(false),
80  inletAreaRatio_(dict.get<scalar>("inletAreaRatio")),
81  inletDischargeCoefficient_(dict.get<scalar>("inletDischargeCoefficient")),
82  timeScale_(dict.lookupOrDefault<scalar>("timeScale", 0)),
83  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
84  UName_(dict.lookupOrDefault<word>("U", "U"))
85 {
86  hasRho_ = dict.readIfPresent("rho", rho_);
87 }
88 
89 
91 (
93  const fvPatch& p,
95  const fvPatchFieldMapper& mapper
96 )
97 :
98  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
99  gamma_(ptf.gamma_),
100  R_(ptf.R_),
101  supplyMassFlowRate_(ptf.supplyMassFlowRate_),
102  supplyTotalTemperature_(ptf.supplyTotalTemperature_),
103  plenumVolume_(ptf.plenumVolume_),
104  plenumDensity_(ptf.plenumDensity_),
105  plenumTemperature_(ptf.plenumTemperature_),
106  rho_(ptf.rho_),
107  hasRho_(ptf.hasRho_),
108  inletAreaRatio_(ptf.inletAreaRatio_),
109  inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
110  timeScale_(ptf.timeScale_),
111  phiName_(ptf.phiName_),
112  UName_(ptf.UName_)
113 {}
114 
115 
117 (
119 )
120 :
121  fixedValueFvPatchScalarField(tppsf),
122  gamma_(tppsf.gamma_),
123  R_(tppsf.R_),
124  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
125  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
126  plenumVolume_(tppsf.plenumVolume_),
127  plenumDensity_(tppsf.plenumDensity_),
128  plenumTemperature_(tppsf.plenumTemperature_),
129  rho_(tppsf.rho_),
130  hasRho_(tppsf.hasRho_),
131  inletAreaRatio_(tppsf.inletAreaRatio_),
132  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
133  timeScale_(tppsf.timeScale_),
134  phiName_(tppsf.phiName_),
135  UName_(tppsf.UName_)
136 {}
137 
138 
140 (
143 )
144 :
145  fixedValueFvPatchScalarField(tppsf, iF),
146  gamma_(tppsf.gamma_),
147  R_(tppsf.R_),
148  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
149  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
150  plenumVolume_(tppsf.plenumVolume_),
151  plenumDensity_(tppsf.plenumDensity_),
152  plenumTemperature_(tppsf.plenumTemperature_),
153  rho_(tppsf.rho_),
154  hasRho_(tppsf.hasRho_),
155  inletAreaRatio_(tppsf.inletAreaRatio_),
156  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
157  timeScale_(tppsf.timeScale_),
158  phiName_(tppsf.phiName_),
159  UName_(tppsf.UName_)
160 {}
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
166 {
167  if (updated())
168  {
169  return;
170  }
171 
172  // Patch properties
173  const fvPatchField<scalar>& p = *this;
174  const fvPatchField<scalar>& p_old =
175  db().lookupObject<volScalarField>
176  (
177  internalField().name()
178  ).oldTime().boundaryField()[patch().index()];
179  const fvPatchField<vector>& U =
180  patch().lookupPatchField<volVectorField, vector>(UName_);
181  const fvsPatchField<scalar>& phi =
182  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
183 
184  // Get the timestep
185  const scalar dt = db().time().deltaTValue();
186 
187  // Check if operating at a new time index and update the old-time properties
188  // if so
189  if (timeIndex_ != db().time().timeIndex())
190  {
191  timeIndex_ = db().time().timeIndex();
192  plenumDensityOld_ = plenumDensity_;
193  plenumTemperatureOld_ = plenumTemperature_;
194  }
195 
196  // Calculate the current mass flow rate
197  scalar massFlowRate(1.0);
198  if (phi.internalField().dimensions() == dimVelocity*dimArea)
199  {
200  if (hasRho_)
201  {
202  massFlowRate = - gSum(rho_*phi);
203  }
204  else
205  {
207  << "The density must be specified when using a volumetric flux."
208  << exit(FatalError);
209  }
210  }
211  else if
212  (
213  phi.internalField().dimensions()
215  )
216  {
217  if (hasRho_)
218  {
220  << "The density must be not specified when using a mass flux."
221  << exit(FatalError);
222  }
223  else
224  {
225  massFlowRate = - gSum(phi);
226  }
227  }
228  else
229  {
231  << "dimensions of phi are not correct"
232  << "\n on patch " << patch().name()
233  << " of field " << internalField().name()
234  << " in file " << internalField().objectPath() << nl
235  << exit(FatalError);
236  }
237 
238  // Calculate the specific heats
239  const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
240 
241  // Calculate the new plenum properties
242  plenumDensity_ =
243  plenumDensityOld_
244  + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
245  plenumTemperature_ =
246  plenumTemperatureOld_
247  + (dt/(plenumDensity_*cv*plenumVolume_))
248  * (
249  supplyMassFlowRate_
250  *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
251  - massFlowRate*R_*plenumTemperature_
252  );
253  const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
254 
255  // Squared velocity magnitude at exit of channels
256  const scalarField U_e(magSqr(U/inletAreaRatio_));
257 
258  // Exit temperature to plenum temperature ratio
259  const scalarField r
260  (
261  1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
262  );
263 
264  // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
265  const scalarField a
266  (
267  (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
268  );
269 
270  // Isentropic exit temperature to plenum temperature ratio
271  const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
272 
273  // Exit pressure to plenum pressure ratio
274  const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
275 
276  // Limit to prevent outflow
277  const scalarField p_new
278  (
279  (1.0 - pos0(phi))*t*plenumPressure + pos0(phi)*max(p, plenumPressure)
280  );
281 
282  // Relaxation fraction
283  const scalar oneByFraction = timeScale_/dt;
284  const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
285 
286  // Set the new value
287  operator==((1.0 - fraction)*p_old + fraction*p_new);
288  fixedValueFvPatchScalarField::updateCoeffs();
289 }
290 
291 
293 {
295  os.writeEntry("gamma", gamma_);
296  os.writeEntry("R", R_);
297  os.writeEntry("supplyMassFlowRate", supplyMassFlowRate_);
298  os.writeEntry("supplyTotalTemperature", supplyTotalTemperature_);
299  os.writeEntry("plenumVolume", plenumVolume_);
300  os.writeEntry("plenumDensity", plenumDensity_);
301  os.writeEntry("plenumTemperature", plenumTemperature_);
302  if (hasRho_)
303  {
304  os.writeEntry("rho", rho_);
305  }
306  os.writeEntry("inletAreaRatio", inletAreaRatio_);
307  os.writeEntry("inletDischargeCoefficient", inletDischargeCoefficient_);
308  os.writeEntryIfDifferent<scalar>("timeScale", 0.0, timeScale_);
309  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
310  os.writeEntryIfDifferent<word>("U", "U", UName_);
311  writeEntry("value", os);
312 }
313 
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 namespace Foam
318 {
320  (
323  );
324 }
325 
326 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:231
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
oldTime
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Definition: createK.H:12
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
plenumPressureFvPatchScalarField.H
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::plenumPressureFvPatchScalarField
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
Definition: plenumPressureFvPatchScalarField.H:193
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::cp
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: MSwindows.C:802
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::plenumPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: plenumPressureFvPatchScalarField.C:292
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::plenumPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: plenumPressureFvPatchScalarField.C:165
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::plenumPressureFvPatchScalarField::plenumPressureFvPatchScalarField
plenumPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: plenumPressureFvPatchScalarField.C:37
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54