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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchScalarField(p, iF),
44  gamma_(1.4),
45  R_(287.04),
46  supplyMassFlowRate_(1.0),
47  supplyTotalTemperature_(300.0),
48  plenumVolume_(1.0),
49  plenumDensity_(1.0),
50  plenumDensityOld_(1.0),
51  plenumTemperature_(300.0),
52  plenumTemperatureOld_(300.0),
53  rho_(1.0),
54  hasRho_(false),
55  inletAreaRatio_(1.0),
56  inletDischargeCoefficient_(1.0),
57  timeScale_(0.0),
58  timeIndex_(-1),
59  phiName_("phi"),
60  UName_("U")
61 {}
62 
63 
65 (
66  const fvPatch& p,
68  const dictionary& dict
69 )
70 :
71  fixedValueFvPatchScalarField(p, iF, dict),
72  gamma_(dict.get<scalar>("gamma")),
73  R_(dict.get<scalar>("R")),
74  supplyMassFlowRate_(dict.get<scalar>("supplyMassFlowRate")),
75  supplyTotalTemperature_(dict.get<scalar>("supplyTotalTemperature")),
76  plenumVolume_(dict.get<scalar>("plenumVolume")),
77  plenumDensity_(dict.get<scalar>("plenumDensity")),
78  plenumTemperature_(dict.get<scalar>("plenumTemperature")),
79  rho_(1.0),
80  hasRho_(false),
81  inletAreaRatio_(dict.get<scalar>("inletAreaRatio")),
82  inletDischargeCoefficient_(dict.get<scalar>("inletDischargeCoefficient")),
83  timeScale_(dict.getOrDefault<scalar>("timeScale", 0)),
84  phiName_(dict.getOrDefault<word>("phi", "phi")),
85  UName_(dict.getOrDefault<word>("U", "U"))
86 {
87  hasRho_ = dict.readIfPresent("rho", rho_);
88 }
89 
90 
92 (
94  const fvPatch& p,
96  const fvPatchFieldMapper& mapper
97 )
98 :
99  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
100  gamma_(ptf.gamma_),
101  R_(ptf.R_),
102  supplyMassFlowRate_(ptf.supplyMassFlowRate_),
103  supplyTotalTemperature_(ptf.supplyTotalTemperature_),
104  plenumVolume_(ptf.plenumVolume_),
105  plenumDensity_(ptf.plenumDensity_),
106  plenumTemperature_(ptf.plenumTemperature_),
107  rho_(ptf.rho_),
108  hasRho_(ptf.hasRho_),
109  inletAreaRatio_(ptf.inletAreaRatio_),
110  inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
111  timeScale_(ptf.timeScale_),
112  phiName_(ptf.phiName_),
113  UName_(ptf.UName_)
114 {}
115 
116 
118 (
120 )
121 :
122  fixedValueFvPatchScalarField(tppsf),
123  gamma_(tppsf.gamma_),
124  R_(tppsf.R_),
125  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
126  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
127  plenumVolume_(tppsf.plenumVolume_),
128  plenumDensity_(tppsf.plenumDensity_),
129  plenumTemperature_(tppsf.plenumTemperature_),
130  rho_(tppsf.rho_),
131  hasRho_(tppsf.hasRho_),
132  inletAreaRatio_(tppsf.inletAreaRatio_),
133  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
134  timeScale_(tppsf.timeScale_),
135  phiName_(tppsf.phiName_),
136  UName_(tppsf.UName_)
137 {}
138 
139 
141 (
144 )
145 :
146  fixedValueFvPatchScalarField(tppsf, iF),
147  gamma_(tppsf.gamma_),
148  R_(tppsf.R_),
149  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
150  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
151  plenumVolume_(tppsf.plenumVolume_),
152  plenumDensity_(tppsf.plenumDensity_),
153  plenumTemperature_(tppsf.plenumTemperature_),
154  rho_(tppsf.rho_),
155  hasRho_(tppsf.hasRho_),
156  inletAreaRatio_(tppsf.inletAreaRatio_),
157  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
158  timeScale_(tppsf.timeScale_),
159  phiName_(tppsf.phiName_),
160  UName_(tppsf.UName_)
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
167 {
168  if (updated())
169  {
170  return;
171  }
172 
173  // Patch properties
174  const fvPatchField<scalar>& p = *this;
175  const fvPatchField<scalar>& p_old =
176  db().lookupObject<volScalarField>
177  (
178  internalField().name()
179  ).oldTime().boundaryField()[patch().index()];
180 
181  const auto& U = patch().lookupPatchField<volVectorField>(UName_);
182  const auto& phi = patch().lookupPatchField<surfaceScalarField>(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() == dimVolume/dimTime)
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 (phi.internalField().dimensions() == dimMass/dimTime)
212  {
213  if (hasRho_)
214  {
216  << "The density must be not specified when using a mass flux."
217  << exit(FatalError);
218  }
219  else
220  {
221  massFlowRate = - gSum(phi);
222  }
223  }
224  else
225  {
227  << "dimensions of phi are not correct"
228  << "\n on patch " << patch().name()
229  << " of field " << internalField().name()
230  << " in file " << internalField().objectPath() << nl
231  << exit(FatalError);
232  }
233 
234  // Calculate the specific heats
235  const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
236 
237  // Calculate the new plenum properties
238  plenumDensity_ =
239  plenumDensityOld_
240  + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
241  plenumTemperature_ =
242  plenumTemperatureOld_
243  + (dt/(plenumDensity_*cv*plenumVolume_))
244  * (
245  supplyMassFlowRate_
246  *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
247  - massFlowRate*R_*plenumTemperature_
248  );
249  const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
250 
251  // Squared velocity magnitude at exit of channels
252  const scalarField U_e(magSqr(U/inletAreaRatio_));
253 
254  // Exit temperature to plenum temperature ratio
255  const scalarField r
256  (
257  1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
258  );
259 
260  // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
261  const scalarField a
262  (
263  (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
264  );
265 
266  // Isentropic exit temperature to plenum temperature ratio
267  const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
268 
269  // Exit pressure to plenum pressure ratio
270  const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
271 
272  // Limit to prevent outflow
273  const scalarField p_new
274  (
275  lerp
276  (
277  t*plenumPressure, // Negative phi
278  max(p, plenumPressure), // Positive phi
279  pos0(phi) // 0-1 selector
280  )
281  );
282 
283  // Relaxation fraction
284  const scalar oneByFraction = timeScale_/dt;
285  const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
287  // Set the new value
288  operator==(lerp(p_old, p_new, fraction));
289  fixedValueFvPatchScalarField::updateCoeffs();
290 }
291 
292 
294 {
296  os.writeEntry("gamma", gamma_);
297  os.writeEntry("R", R_);
298  os.writeEntry("supplyMassFlowRate", supplyMassFlowRate_);
299  os.writeEntry("supplyTotalTemperature", supplyTotalTemperature_);
300  os.writeEntry("plenumVolume", plenumVolume_);
301  os.writeEntry("plenumDensity", plenumDensity_);
302  os.writeEntry("plenumTemperature", plenumTemperature_);
303  if (hasRho_)
304  {
305  os.writeEntry("rho", rho_);
306  }
307  os.writeEntry("inletAreaRatio", inletAreaRatio_);
308  os.writeEntry("inletDischargeCoefficient", inletDischargeCoefficient_);
309  os.writeEntryIfDifferent<scalar>("timeScale", 0.0, timeScale_);
310  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
311  os.writeEntryIfDifferent<word>("U", "U", UName_);
313 }
314 
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 namespace Foam
319 {
321  (
323  plenumPressureFvPatchScalarField
324  );
325 }
326 
327 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
dictionary dict
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
Info<< "Creating field kinetic energy K\"<< 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
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: POSIX.C:1063
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Type gSum(const FieldField< Field, Type > &f)
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
plenumPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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:336
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition: IOobject.C:456
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
dimensionedScalar pos0(const dimensionedScalar &ds)
OBJstream os(runTime.globalPath()/outputName)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
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))
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24