flowRateInletVelocityFvPatchVectorField.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) 2011-2024 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 "volFields.H"
29 #include "one.H"
30 #include "patchPatchDist.H"
31 #include "wallPolyPatch.H"
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 void Foam::flowRateInletVelocityFvPatchVectorField::setWallDist()
36 {
37  if (profile_.valid())
38  {
39  const labelHashSet otherPatchIDs
40  (
41  patch().patch().boundaryMesh().findIndices<wallPolyPatch>()
42  );
43 
44  const patchPatchDist pwd(patch().patch(), otherPatchIDs);
45 
46  y_ = pwd/gMax(pwd);
47  }
48 
49  area_ = gSum(patch().magSf());
50 }
51 
52 
54 Foam::flowRateInletVelocityFvPatchVectorField::profile()
55 {
56  if (profile_.valid())
57  {
58  return profile_->value(y_);
59  }
60  else
61  {
62  return tmp<scalarField>(new scalarField(size(), scalar(1)));
63  }
64 }
65 
66 
67 template<class ScaleType, class AlphaType, class RhoType>
68 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
69 (
70  const ScaleType& scale,
71  const AlphaType& alpha,
72  const RhoType& rho
73 )
74 {
75  const scalarField profile(this->profile());
76 
77  const scalar avgU =
78  scale
79  *flowRate_->value(db().time().value())
80  /gSum(alpha*rho*profile*patch().magSf());
81 
82  operator==(- avgU*profile*patch().nf());
83 }
84 
85 
86 template<class AlphaType>
87 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
88 (
89  const AlphaType& alpha
90 )
91 {
92  if (meanVelocity_)
93  {
94  updateValues(area_, alpha, one());
95  }
96  else if (volumetric_ || rhoName_ == "none")
97  {
98  updateValues(one(), alpha, one());
99  }
100  else
101  {
102  // Mass flow-rate
103  if (db().foundObject<volScalarField>(rhoName_))
104  {
105  const fvPatchField<scalar>& rhop =
106  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
107 
108  updateValues(one(), alpha, rhop);
109  }
110  else
111  {
112  // Use constant density
113  if (rhoInlet_ < 0)
114  {
116  << "Did not find registered density field " << rhoName_
117  << " and no constant density 'rhoInlet' specified"
118  << exit(FatalError);
119  }
120 
121  updateValues(one(), alpha, rhoInlet_);
122  }
123  }
124 }
125 
126 
127 bool Foam::flowRateInletVelocityFvPatchVectorField::canEvaluate()
128 {
129  return
131  || !patch().boundaryMesh().mesh().time().processorCase();
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
136 
139 (
140  const fvPatch& p,
142  const dictionary& dict
143 )
144 :
145  fixedValueFvPatchField<vector>(p, iF, dict, false),
146  flowRate_(),
147  meanVelocity_(),
148  volumetric_(),
149  profile_(),
150  rhoName_("rho"),
151  rhoInlet_(dict.lookupOrDefault<scalar>("rhoInlet", dimDensity, -vGreat)),
152  alphaName_(dict.lookupOrDefault<word>("alpha", word::null)),
153  y_(),
154  area_(NaN)
155 {
156  if (dict.found("meanVelocity"))
157  {
158  meanVelocity_ = true;
159  volumetric_ = false;
160  flowRate_ =
162  (
163  "meanVelocity",
164  db().time().userUnits(),
165  dimVelocity,
166  dict
167  );
168  }
169  else if (dict.found("volumetricFlowRate"))
170  {
171  meanVelocity_ = false;
172  volumetric_ = true;
173  flowRate_ =
175  (
176  "volumetricFlowRate",
177  db().time().userUnits(),
179  dict
180  );
181  }
182  else if (dict.found("massFlowRate"))
183  {
184  meanVelocity_ = false;
185  volumetric_ = false;
186  flowRate_ =
188  (
189  "massFlowRate",
190  db().time().userUnits(),
191  dimMassFlux,
192  dict
193  );
194  rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
195  }
196  else
197  {
199  << "Please supply 'meanVelocity', 'volumetricFlowRate' or"
200  << " 'massFlowRate' and 'rho'" << exit(FatalIOError);
201  }
202 
203  if (dict.found("profile"))
204  {
205  profile_ = Function1<scalar>::New("profile", dimless, dimless, dict);
206  }
207 
208  if (canEvaluate())
209  {
210  setWallDist();
211  }
212 
213  if (!canEvaluate() || dict.found("value"))
214  {
216  (
217  vectorField("value", iF.dimensions(), dict, p.size())
218  );
219  }
220  else
221  {
223  }
224 }
225 
226 
229 (
231  const fvPatch& p,
233  const fieldMapper& mapper
234 )
235 :
236  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
237  flowRate_(ptf.flowRate_, false),
238  meanVelocity_(ptf.meanVelocity_),
239  volumetric_(ptf.volumetric_),
240  profile_(ptf.profile_, false),
241  rhoName_(ptf.rhoName_),
242  rhoInlet_(ptf.rhoInlet_),
243  alphaName_(ptf.alphaName_),
244  y_
245  (
246  profile_.valid() && canEvaluate()
247  ? mapper(ptf.y_)
248  : tmp<scalarField>(new scalarField())
249  ),
250  area_(NaN)
251 {}
252 
253 
256 (
259 )
260 :
261  fixedValueFvPatchField<vector>(ptf, iF),
262  flowRate_(ptf.flowRate_, false),
263  meanVelocity_(ptf.meanVelocity_),
264  volumetric_(ptf.volumetric_),
265  profile_(ptf.profile_, false),
266  rhoName_(ptf.rhoName_),
267  rhoInlet_(ptf.rhoInlet_),
268  alphaName_(ptf.alphaName_),
269  y_(ptf.y_),
270  area_(ptf.area_)
271 {}
272 
273 
274 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
275 
277 (
278  const fvPatchVectorField& ptf,
279  const fieldMapper& mapper
280 )
281 {
282  fixedValueFvPatchVectorField::map(ptf, mapper);
283 
285  refCast<const flowRateInletVelocityFvPatchVectorField>(ptf);
286 
287  if (profile_.valid() && canEvaluate())
288  {
289  mapper(y_, tiptf.y_);
290  }
291 }
292 
293 
295 (
296  const fvPatchVectorField& ptf
297 )
298 {
299  fixedValueFvPatchVectorField::reset(ptf);
300 
302  refCast<const flowRateInletVelocityFvPatchVectorField>(ptf);
303 
304  if (profile_.valid() && canEvaluate())
305  {
306  y_.reset(tiptf.y_);
307  }
308 }
309 
310 
312 {
313  if (updated())
314  {
315  return;
316  }
317 
318  if (!canEvaluate())
319  {
321  << "Cannot evaluate flow rate on a non-parallel processor case"
322  << exit(FatalError);
323  }
324 
325  if (alphaName_ != word::null)
326  {
327  const fvPatchField<scalar>& alphap =
328  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
329 
330  updateValues(alphap);
331  }
332  else
333  {
334  updateValues(one());
335  }
336 
337  fixedValueFvPatchVectorField::updateCoeffs();
338 }
339 
340 
342 {
344  writeEntry(os, db().time().userUnits(), unitAny, flowRate_());
345  if (profile_.valid())
346  {
347  writeEntry(os, profile_());
348  }
349  if (!volumetric_)
350  {
351  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
352  writeEntryIfDifferent<scalar>(os, "rhoInlet", -vGreat, rhoInlet_);
353  }
354  writeEntryIfDifferent<word>(os, "alpha", word::null, alphaName_);
355  writeEntry(os, "value", *this);
356 }
357 
358 
359 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
360 
361 namespace Foam
362 {
364  (
367  );
368 }
369 
370 
371 // ************************************************************************* //
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 dimensionSet & dimensions() const
Return dimensions.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Velocity inlet boundary condition creating a velocity field with optionally specified profile normal ...
flowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:51
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
bool valid(const PtrList< ModelType > &l)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
const unitConversion unitAny
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimMassFlux
const dimensionSet dimVolumetricFlux
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
const dimensionSet dimVelocity
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMax(const FieldField< Field, Type > &f)
dictionary dict
volScalarField & p