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-------------------------------------------------------------------------------
11License
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
31#include "fvPatchFieldMapper.H"
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 const fvPatchField<vector>& U =
181 patch().lookupPatchField<volVectorField, vector>(UName_);
183 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
184
185 // Get the timestep
186 const scalar dt = db().time().deltaTValue();
187
188 // Check if operating at a new time index and update the old-time properties
189 // if so
190 if (timeIndex_ != db().time().timeIndex())
191 {
192 timeIndex_ = db().time().timeIndex();
193 plenumDensityOld_ = plenumDensity_;
194 plenumTemperatureOld_ = plenumTemperature_;
195 }
196
197 // Calculate the current mass flow rate
198 scalar massFlowRate(1.0);
200 {
201 if (hasRho_)
202 {
203 massFlowRate = - gSum(rho_*phi);
204 }
205 else
206 {
208 << "The density must be specified when using a volumetric flux."
209 << exit(FatalError);
210 }
211 }
212 else if
213 (
216 )
217 {
218 if (hasRho_)
219 {
221 << "The density must be not specified when using a mass flux."
222 << exit(FatalError);
223 }
224 else
225 {
226 massFlowRate = - gSum(phi);
227 }
228 }
229 else
230 {
232 << "dimensions of phi are not correct"
233 << "\n on patch " << patch().name()
234 << " of field " << internalField().name()
235 << " in file " << internalField().objectPath() << nl
236 << exit(FatalError);
237 }
238
239 // Calculate the specific heats
240 const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
241
242 // Calculate the new plenum properties
243 plenumDensity_ =
244 plenumDensityOld_
245 + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
246 plenumTemperature_ =
247 plenumTemperatureOld_
248 + (dt/(plenumDensity_*cv*plenumVolume_))
249 * (
250 supplyMassFlowRate_
251 *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
252 - massFlowRate*R_*plenumTemperature_
253 );
254 const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
255
256 // Squared velocity magnitude at exit of channels
257 const scalarField U_e(magSqr(U/inletAreaRatio_));
258
259 // Exit temperature to plenum temperature ratio
260 const scalarField r
261 (
262 1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
263 );
264
265 // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
266 const scalarField a
267 (
268 (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
269 );
270
271 // Isentropic exit temperature to plenum temperature ratio
272 const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
273
274 // Exit pressure to plenum pressure ratio
275 const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
276
277 // Limit to prevent outflow
278 const scalarField p_new
279 (
280 (1.0 - pos0(phi))*t*plenumPressure + pos0(phi)*max(p, plenumPressure)
281 );
282
283 // Relaxation fraction
284 const scalar oneByFraction = timeScale_/dt;
285 const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
286
287 // Set the new value
288 operator==((1.0 - fraction)*p_old + fraction*p_new);
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_);
312 writeEntry("value", os);
313}
314
315
316// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317
318namespace Foam
319{
321 (
324 );
325}
326
327// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:506
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
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()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
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))
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensionedScalar sqrt(const dimensionedScalar &ds)
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
const volScalarField & cp
Foam::surfaceFields.