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 -------------------------------------------------------------------------------
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 
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "TableFile.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 const Foam::Enum
38 <
40 >
42 ({
43  { fanFlowDirection::ffdIn, "in" },
44  { fanFlowDirection::ffdOut, "out" },
45 });
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const fvPatch& p,
53  const DimensionedField<scalar, volMesh>& iF
54 )
55 :
56  totalPressureFvPatchScalarField(p, iF),
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 (
119  const fanPressureFvPatchScalarField& rhs
120 )
121 :
122  totalPressureFvPatchScalarField(rhs),
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& phip = patch().lookupPatchField<surfaceScalarField>(phiName());
157 
158  const int dir = 2*direction_ - 1;
159 
160  // Average volumetric flow rate
161  scalar volFlowRate = 0;
162 
163  if (phip.internalField().dimensions() == dimVolume/dimTime)
164  {
165  volFlowRate = dir*gSum(phip);
166  }
167  else if (phip.internalField().dimensions() == dimMass/dimTime)
168  {
169  const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName());
170  volFlowRate = dir*gSum(phip/rhop);
171  }
172  else
173  {
175  << "dimensions of phi are not correct\n"
176  << " on patch " << patch().name()
177  << " of field " << internalField().name()
178  << " in file " << internalField().objectPath() << nl
179  << exit(FatalError);
180  }
181 
182  // The non-dimensional parameters
183 
184  scalar rpm(0);
185  scalar meanDiam(0);
186 
187  if (nonDimensional_)
188  {
189  rpm = rpm_->value(this->db().time().timeOutputValue());
190  meanDiam = dm_->value(this->db().time().timeOutputValue());
191 
192  // Create an non-dimensional flow rate
193  volFlowRate =
194  120.0*volFlowRate
195  / stabilise
196  (
197  pow3(constant::mathematical::pi * meanDiam) * rpm,
198  VSMALL
199  );
200  }
201 
202  // Pressure drop for this flow rate
203  scalar pdFan = fanCurve_->value(max(volFlowRate, scalar(0)));
204 
205  if (nonDimensional_)
206  {
207  // Convert the non-dimensional deltap from curve into deltaP
208  pdFan =
209  (
211  * sqr(rpm * meanDiam) / 1800.0
212  );
213  }
214 
216  (
217  p0() - dir*pdFan,
218  patch().lookupPatchField<volVectorField>(UName())
219  );
220 }
221 
222 
224 {
226  fanCurve_->writeData(os);
227  os.writeEntry("direction", fanFlowDirectionNames_[direction_]);
228 
229  if (nonDimensional_)
230  {
231  os.writeEntry("nonDimensional", "true");
232  rpm_->writeData(os);
233  dm_->writeData(os);
234  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
243  (
245  fanPressureFvPatchScalarField
246  );
247 };
248 
249 
250 // ************************************************************************* //
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
Templated table container function where data are read from file.
Definition: TableFile.H:74
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
const word UName(propsDict.getOrDefault< word >("U", "U"))
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
This boundary condition can be applied to assign either a pressure inlet or outlet total pressure con...
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
static const Enum< fanFlowDirection > fanFlowDirectionNames_
Fan flow direction names.
A FieldMapper for finite-volume patch fields.
constexpr scalar pi(M_PI)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar pow3(const dimensionedScalar &ds)
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
dimensionedScalar pow4(const dimensionedScalar &ds)
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
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const volScalarField & p0
Definition: EEqn.H:36
virtual void write(Ostream &os) const
Write.
Namespace for OpenFOAM.
fanPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
This boundary condition provides a total pressure condition. Four variants are possible: ...