Loading...
Searching...
No Matches
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,2024 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 plenumDensityOld_(1.0),
79 plenumTemperature_(dict.get<scalar>("plenumTemperature")),
80 plenumTemperatureOld_(300.0),
81 rho_(1.0),
82 hasRho_(false),
83 inletAreaRatio_(dict.get<scalar>("inletAreaRatio")),
84 inletDischargeCoefficient_(dict.get<scalar>("inletDischargeCoefficient")),
85 timeScale_(dict.getOrDefault<scalar>("timeScale", 0)),
86 timeIndex_(-1),
87 phiName_(dict.getOrDefault<word>("phi", "phi")),
88 UName_(dict.getOrDefault<word>("U", "U"))
89{
90 hasRho_ = dict.readIfPresent("rho", rho_);
91}
92
93
95(
97 const fvPatch& p,
99 const fvPatchFieldMapper& mapper
100)
101:
102 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
103 gamma_(ptf.gamma_),
104 R_(ptf.R_),
105 supplyMassFlowRate_(ptf.supplyMassFlowRate_),
106 supplyTotalTemperature_(ptf.supplyTotalTemperature_),
107 plenumVolume_(ptf.plenumVolume_),
108 plenumDensity_(ptf.plenumDensity_),
109 plenumDensityOld_(ptf.plenumDensityOld_),
110 plenumTemperature_(ptf.plenumTemperature_),
111 plenumTemperatureOld_(ptf.plenumTemperatureOld_),
112 rho_(ptf.rho_),
113 hasRho_(ptf.hasRho_),
114 inletAreaRatio_(ptf.inletAreaRatio_),
115 inletDischargeCoefficient_(ptf.inletDischargeCoefficient_),
116 timeScale_(ptf.timeScale_),
117 timeIndex_(ptf.timeIndex_),
118 phiName_(ptf.phiName_),
119 UName_(ptf.UName_)
120{}
121
122
124(
126)
127:
128 fixedValueFvPatchScalarField(tppsf),
129 gamma_(tppsf.gamma_),
130 R_(tppsf.R_),
131 supplyMassFlowRate_(tppsf.supplyMassFlowRate_),
132 supplyTotalTemperature_(tppsf.supplyTotalTemperature_),
133 plenumVolume_(tppsf.plenumVolume_),
134 plenumDensity_(tppsf.plenumDensity_),
135 plenumDensityOld_(tppsf.plenumDensityOld_),
136 plenumTemperature_(tppsf.plenumTemperature_),
137 plenumTemperatureOld_(tppsf.plenumTemperatureOld_),
138 rho_(tppsf.rho_),
139 hasRho_(tppsf.hasRho_),
140 inletAreaRatio_(tppsf.inletAreaRatio_),
141 inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
142 timeScale_(tppsf.timeScale_),
143 timeIndex_(tppsf.timeIndex_),
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 plenumDensityOld_(tppsf.plenumDensityOld_),
163 plenumTemperature_(tppsf.plenumTemperature_),
164 plenumTemperatureOld_(tppsf.plenumTemperatureOld_),
165 rho_(tppsf.rho_),
166 hasRho_(tppsf.hasRho_),
167 inletAreaRatio_(tppsf.inletAreaRatio_),
168 inletDischargeCoefficient_(tppsf.inletDischargeCoefficient_),
169 timeScale_(tppsf.timeScale_),
170 timeIndex_(tppsf.timeIndex_),
171 phiName_(tppsf.phiName_),
172 UName_(tppsf.UName_)
173{}
174
175
176// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177
179{
180 if (updated())
181 {
182 return;
183 }
184
185 // Patch properties
186 const fvPatchField<scalar>& p = *this;
187 const fvPatchField<scalar>& p_old =
188 db().lookupObject<volScalarField>
189 (
190 internalField().name()
191 ).oldTime().boundaryField()[patch().index()];
192
193 const auto& U = patch().lookupPatchField<volVectorField>(UName_);
194 const auto& phi = patch().lookupPatchField<surfaceScalarField>(phiName_);
195
196 // Get the timestep
197 const scalar dt = db().time().deltaTValue();
198
199 // Check if operating at a new time index and update the old-time properties
200 // if so
201 if (timeIndex_ != db().time().timeIndex())
202 {
203 timeIndex_ = db().time().timeIndex();
204 plenumDensityOld_ = plenumDensity_;
205 plenumTemperatureOld_ = plenumTemperature_;
206 }
207
208 // Calculate the current mass flow rate
209 scalar massFlowRate(1.0);
210 if (phi.internalField().dimensions() == dimVolume/dimTime)
211 {
212 if (hasRho_)
213 {
214 massFlowRate = - rho_*gSum(phi);
215 }
216 else
217 {
219 << "The density must be specified when using a volumetric flux."
220 << exit(FatalError);
221 }
222 }
223 else if (phi.internalField().dimensions() == dimMass/dimTime)
224 {
225 if (hasRho_)
226 {
228 << "The density must be not specified when using a mass flux."
229 << exit(FatalError);
230 }
231 else
232 {
233 massFlowRate = - gSum(phi);
234 }
235 }
236 else
237 {
239 << "dimensions of phi are not correct"
240 << "\n on patch " << patch().name()
241 << " of field " << internalField().name()
242 << " in file " << internalField().objectPath() << nl
243 << exit(FatalError);
244 }
245
246 // Calculate the specific heats
247 const scalar cv = R_/(gamma_ - 1), cp = R_*gamma_/(gamma_ - 1);
248
249 // Calculate the new plenum properties
250 plenumDensity_ =
251 plenumDensityOld_
252 + (dt/plenumVolume_)*(supplyMassFlowRate_ - massFlowRate);
253 plenumTemperature_ =
254 plenumTemperatureOld_
255 + (dt/(plenumDensity_*cv*plenumVolume_))
256 * (
257 supplyMassFlowRate_
258 *(cp*supplyTotalTemperature_ - cv*plenumTemperature_)
259 - massFlowRate*R_*plenumTemperature_
260 );
261 const scalar plenumPressure = plenumDensity_*R_*plenumTemperature_;
262
263 // Squared velocity magnitude at exit of channels
264 const scalarField U_e(magSqr(U/inletAreaRatio_));
265
266 // Exit temperature to plenum temperature ratio
267 const scalarField r
268 (
269 1.0 - (gamma_ - 1.0)*U_e/(2.0*gamma_*R_*plenumTemperature_)
270 );
271
272 // Quadratic coefficient (others not needed as b = +1.0 and c = -1.0)
273 const scalarField a
274 (
275 (1.0 - r)/(r*r*inletDischargeCoefficient_*inletDischargeCoefficient_)
276 );
277
278 // Isentropic exit temperature to plenum temperature ratio
279 const scalarField s(2.0/(1.0 + sqrt(1.0 + 4.0*a)));
280
281 // Exit pressure to plenum pressure ratio
282 const scalarField t(pow(s, gamma_/(gamma_ - 1.0)));
283
284 // Limit to prevent outflow
285 const scalarField p_new
286 (
287 lerp
288 (
289 t*plenumPressure, // Negative phi
290 max(p, plenumPressure), // Positive phi
291 pos0(phi) // 0-1 selector
292 )
293 );
294
295 // Relaxation fraction
296 const scalar oneByFraction = timeScale_/dt;
297 const scalar fraction = oneByFraction < 1.0 ? 1.0 : 1.0/oneByFraction;
299 // Set the new value
300 operator==(lerp(p_old, p_new, fraction));
301 fixedValueFvPatchScalarField::updateCoeffs();
302}
303
304
306{
308 os.writeEntry("gamma", gamma_);
309 os.writeEntry("R", R_);
310 os.writeEntry("supplyMassFlowRate", supplyMassFlowRate_);
311 os.writeEntry("supplyTotalTemperature", supplyTotalTemperature_);
312 os.writeEntry("plenumVolume", plenumVolume_);
313 os.writeEntry("plenumDensity", plenumDensity_);
314 os.writeEntry("plenumTemperature", plenumTemperature_);
315 if (hasRho_)
316 {
317 os.writeEntry("rho", rho_);
318 }
319 os.writeEntry("inletAreaRatio", inletAreaRatio_);
320 os.writeEntry("inletDischargeCoefficient", inletDischargeCoefficient_);
321 os.writeEntryIfDifferent<scalar>("timeScale", 0.0, timeScale_);
322 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
323 os.writeEntryIfDifferent<word>("U", "U", UName_);
326
327
328// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329
330namespace Foam
331{
333 (
335 plenumPressureFvPatchScalarField
336 );
337}
338
339// ************************************************************************* //
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...
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
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:59
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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...
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a plenum pressure inlet condition. This condition creates a zero-dim...
plenumPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition word.H:66
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:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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))
const std::string patch
OpenFOAM patch number as a std::string.
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:40
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
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:1065
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
dictionary dict
const volScalarField & cp
Foam::surfaceFields.