fanPressureFvPatchScalarField.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2021 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 "volFields.H"
32#include "surfaceFields.H"
33#include "TableFile.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37const Foam::Enum
38<
40>
42({
43 { fanFlowDirection::ffdIn, "in" },
44 { fanFlowDirection::ffdOut, "out" },
45});
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const fvPatch& p,
54)
55:
57 fanCurve_(nullptr),
58 direction_(ffdOut),
59 nonDimensional_(false),
60 rpm_(nullptr),
61 dm_(nullptr)
62{}
63
64
66(
68 const fvPatch& p,
70 const fvPatchFieldMapper& mapper
71)
72:
73 totalPressureFvPatchScalarField(rhs, p, iF, mapper),
74 fanCurve_(rhs.fanCurve_.clone()),
75 direction_(rhs.direction_),
76 nonDimensional_(rhs.nonDimensional_),
77 rpm_(rhs.rpm_.clone()),
78 dm_(rhs.dm_.clone())
79{}
80
81
83(
84 const fvPatch& p,
86 const dictionary& dict
87)
88:
90 fanCurve_(nullptr),
91 direction_(fanFlowDirectionNames_.get("direction", dict)),
92 nonDimensional_(dict.getOrDefault("nonDimensional", false)),
93 rpm_(nullptr),
94 dm_(nullptr)
95{
96 // Backwards compatibility
97 if (dict.found("file"))
98 {
99 fanCurve_.reset
100 (
101 new Function1Types::TableFile<scalar>("fanCurve", dict, &this->db())
102 );
103 }
104 else
105 {
106 fanCurve_.reset(Function1<scalar>::New("fanCurve", dict, &this->db()));
107 }
108
109 if (nonDimensional_)
110 {
111 rpm_.reset(Function1<scalar>::New("rpm", dict, &this->db()));
112 dm_.reset(Function1<scalar>::New("dm", dict, &this->db()));
113 }
114}
115
116
118(
120)
121:
123 fanCurve_(rhs.fanCurve_.clone()),
124 direction_(rhs.direction_),
125 nonDimensional_(rhs.nonDimensional_),
126 rpm_(rhs.rpm_.clone()),
127 dm_(rhs.dm_.clone())
128{}
129
130
132(
135)
136:
138 fanCurve_(rhs.fanCurve_.clone()),
139 direction_(rhs.direction_),
140 nonDimensional_(rhs.nonDimensional_),
141 rpm_(rhs.rpm_.clone()),
142 dm_(rhs.dm_.clone())
143{}
144
145
146// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147
149{
150 if (updated())
151 {
152 return;
153 }
154
155 // Retrieve flux field
156 const auto& phi = db().lookupObject<surfaceScalarField>(phiName());
157
158 const auto& phip = patch().patchField<surfaceScalarField, scalar>(phi);
159
160 const int dir = 2*direction_ - 1;
161
162 // Average volumetric flow rate
163 scalar volFlowRate = 0;
164
166 {
167 volFlowRate = dir*gSum(phip);
168 }
170 {
171 const scalarField& rhop =
172 patch().lookupPatchField<volScalarField, scalar>(rhoName());
173 volFlowRate = dir*gSum(phip/rhop);
174 }
175 else
176 {
178 << "dimensions of phi are not correct\n"
179 << " on patch " << patch().name()
180 << " of field " << internalField().name()
181 << " in file " << internalField().objectPath() << nl
182 << exit(FatalError);
183 }
184
185 // The non-dimensional parameters
186
187 scalar rpm(0);
188 scalar meanDiam(0);
189
190 if (nonDimensional_)
191 {
192 rpm = rpm_->value(this->db().time().timeOutputValue());
193 meanDiam = dm_->value(this->db().time().timeOutputValue());
194
195 // Create an non-dimensional flow rate
196 volFlowRate =
197 120.0*volFlowRate
198 / stabilise
199 (
200 pow3(constant::mathematical::pi * meanDiam) * rpm,
201 VSMALL
202 );
203 }
204
205 // Pressure drop for this flow rate
206 scalar pdFan = fanCurve_->value(max(volFlowRate, scalar(0)));
207
208 if (nonDimensional_)
209 {
210 // Convert the non-dimensional deltap from curve into deltaP
211 pdFan =
212 (
214 * sqr(rpm * meanDiam) / 1800.0
215 );
216 }
217
219 (
220 p0() - dir*pdFan,
221 patch().lookupPatchField<volVectorField, vector>(UName())
222 );
223}
224
225
227{
229 fanCurve_->writeData(os);
230 os.writeEntry("direction", fanFlowDirectionNames_[direction_]);
231
232 if (nonDimensional_)
233 {
234 os.writeEntry("nonDimensional", "true");
235 rpm_->writeData(os);
236 dm_->writeData(os);
237 }
238}
239
240
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243namespace Foam
244{
246 (
249 );
250};
251
252
253// ************************************************************************* //
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.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Templated table container function where data are read from file.
Definition: TableFile.H:81
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
This boundary condition can be applied to assign either a pressure inlet or outlet total pressure con...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static const Enum< fanFlowDirection > fanFlowDirectionNames_
Fan flow direction names.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a total pressure condition. Four variants are possible:
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
const volScalarField & p0
Definition: EEqn.H:36
#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
constexpr scalar pi(M_PI)
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)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
Foam::surfaceFields.
const word UName(propsDict.getOrDefault< word >("U", "U"))