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-2023 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  const dictionary& dict
39 )
40 :
41  fixedValueFvPatchScalarField(p, iF, dict),
42  gamma_(dict.lookup<scalar>("gamma")),
43  R_(dict.lookup<scalar>("R")),
44  supplyMassFlowRate_(dict.lookup<scalar>("supplyMassFlowRate")),
45  supplyTotalTemperature_
46  (
47  dict.lookup<scalar>("supplyTotalTemperature")
48  ),
49  plenumVolume_(dict.lookup<scalar>("plenumVolume")),
50  plenumDensity_(dict.lookup<scalar>("plenumDensity")),
51  plenumTemperature_(dict.lookup<scalar>("plenumTemperature")),
52  rho_(1.0),
53  hasRho_(false),
54  inletAreaRatio_(dict.lookup<scalar>("inletAreaRatio")),
55  inletDischargeCoefficient_
56  (
57  dict.lookup<scalar>("inletDischargeCoefficient")
58  ),
59  timeScale_(dict.lookupOrDefault<scalar>("timeScale", 0.0)),
60  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
61  UName_(dict.lookupOrDefault<word>("U", "U"))
62 {
63  if (dict.found("rho"))
64  {
65  rho_ = dict.lookup<scalar>("rho");
66  hasRho_ = true;
67  }
68 }
69 
70 
72 (
74  const fvPatch& p,
76  const fvPatchFieldMapper& mapper
77 )
78 :
79  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
80  gamma_(ptf.gamma_),
81  R_(ptf.R_),
82  supplyMassFlowRate_(ptf.supplyMassFlowRate_),
83  supplyTotalTemperature_(ptf.supplyTotalTemperature_),
84  plenumVolume_(ptf.plenumVolume_),
85  plenumDensity_(ptf.plenumDensity_),
86  plenumTemperature_(ptf.plenumTemperature_),
87  rho_(ptf.rho_),
88  hasRho_(ptf.hasRho_),
89  inletAreaRatio_(ptf.inletAreaRatio_),
90  inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
91  timeScale_(ptf.timeScale_),
92  phiName_(ptf.phiName_),
93  UName_(ptf.UName_)
94 {}
95 
96 
98 (
101 )
102 :
103  fixedValueFvPatchScalarField(tppsf, iF),
104  gamma_(tppsf.gamma_),
105  R_(tppsf.R_),
106  supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
107  supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
108  plenumVolume_(tppsf.plenumVolume_),
109  plenumDensity_(tppsf.plenumDensity_),
110  plenumTemperature_(tppsf.plenumTemperature_),
111  rho_(tppsf.rho_),
112  hasRho_(tppsf.hasRho_),
113  inletAreaRatio_(tppsf.inletAreaRatio_),
114  inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
115  timeScale_(tppsf.timeScale_),
116  phiName_(tppsf.phiName_),
117  UName_(tppsf.UName_)
118 {}
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 {
125  if (updated())
126  {
127  return;
128  }
129 
130  // Patch properties
131  const fvPatchField<scalar>& p = *this;
132  const fvPatchField<scalar>& p_old =
133  db().lookupObject<volScalarField>
134  (
135  internalField().name()
136  ).oldTime().boundaryField()[patch().index()];
137  const fvPatchField<vector>& U =
138  patch().lookupPatchField<volVectorField, vector>(UName_);
139  const fvsPatchField<scalar>& phi =
140  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
141 
142  // Get the timestep
143  const scalar dt = db().time().deltaTValue();
144 
145  // Check if operating at a new time index and update the old-time properties
146  // if so
147  if (timeIndex_ != db().time().timeIndex())
148  {
149  timeIndex_ = db().time().timeIndex();
150  plenumDensityOld_ = plenumDensity_;
151  plenumTemperatureOld_ = plenumTemperature_;
152  }
153 
154  // Calculate the current mass flow rate
155  scalar massFlowRate(1.0);
156  if (phi.internalField().dimensions() == dimFlux)
157  {
158  if (hasRho_)
159  {
160  massFlowRate = - gSum(rho_*phi);
161  }
162  else
163  {
165  << "The density must be specified when using a volumetric flux."
166  << exit(FatalError);
167  }
168  }
169  else if
170  (
171  phi.internalField().dimensions()
172  == dimMassFlux
173  )
174  {
175  if (hasRho_)
176  {
178  << "The density must be not specified when using a mass flux."
179  << exit(FatalError);
180  }
181  else
182  {
183  massFlowRate = - gSum(phi);
184  }
185  }
186  else
187  {
189  << "dimensions of phi are not correct"
190  << "\n on patch " << patch().name()
191  << " of field " << internalField().name()
192  << " in file " << internalField().objectPath() << nl
193  << exit(FatalError);
194  }
195 
196  // Calcaulate the specific heats
197  const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
198 
199  // Calculate the new plenum properties
200  plenumDensity_ =
201  plenumDensityOld_
202  + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
203  plenumTemperature_ =
204  plenumTemperatureOld_
205  + (dt/(plenumDensity_*cv*plenumVolume_))
206  * (
207  supplyMassFlowRate_
208  *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
209  - massFlowRate*R_*plenumTemperature_
210  );
211  const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
212 
213  // Squared velocity magnitude at exit of channels
214  const scalarField U_e(magSqr(U/inletAreaRatio_));
215 
216  // Exit temperature to plenum temperature ratio
217  const scalarField r
218  (
219  1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
220  );
221 
222  // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
223  const scalarField a
224  (
225  (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
226  );
227 
228  // Isentropic exit temperature to plenum temperature ratio
229  const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
230 
231  // Exit pressure to plenum pressure ratio
232  const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
233 
234  // Limit to prevent outflow
235  const scalarField p_new
236  (
237  neg(phi)*t*plenumPressure + pos0(phi)*max(p, plenumPressure)
238  );
239 
240  // Relaxation fraction
241  const scalar oneByFraction = timeScale_/dt;
242  const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
243 
244  // Set the new value
245  operator==((1.0 - fraction)*p_old + fraction*p_new);
246  fixedValueFvPatchScalarField::updateCoeffs();
247 }
248 
249 
251 {
253  writeEntry(os, "gamma", gamma_);
254  writeEntry(os, "R", R_);
255  writeEntry(os, "supplyMassFlowRate", supplyMassFlowRate_);
256  writeEntry(os, "supplyTotalTemperature", supplyTotalTemperature_);
257  writeEntry(os, "plenumVolume", plenumVolume_);
258  writeEntry(os, "plenumDensity", plenumDensity_);
259  writeEntry(os, "plenumTemperature", plenumTemperature_);
260  if (hasRho_)
261  {
262  writeEntry(os, "rho", rho_);
263  }
264  writeEntry(os, "inletAreaRatio", inletAreaRatio_);
265  writeEntry
266  (
267  os,
268  "inletDischargeCoefficient",
269  inletDischargeCoefficient_
270  );
271  writeEntryIfDifferent<scalar>(os, "timeScale", 0.0, timeScale_);
272  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
273  writeEntryIfDifferent<word>(os, "U", "U", UName_);
274  writeEntry(os, "value", *this);
275 }
276 
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 
280 namespace Foam
281 {
283  (
286  );
287 }
288 
289 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
plenumPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar neg(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
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const dimensionSet dimFlux
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
volScalarField & p
Foam::surfaceFields.