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-2020 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 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38 )
39 :
41  flowRate_(),
42  volumetric_(false),
43  rhoName_("rho"),
44  rhoInlet_(0.0),
45  alphaName_(word::null),
46  extrapolateProfile_(false)
47 {}
48 
49 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
59  rhoInlet_(dict.lookupOrDefault<scalar>("rhoInlet", -vGreat)),
60  alphaName_(dict.lookupOrDefault<word>("alpha", word::null)),
61  extrapolateProfile_
62  (
63  dict.lookupOrDefault<Switch>("extrapolateProfile", false)
64  )
65 {
66  if (dict.found("volumetricFlowRate"))
67  {
68  volumetric_ = true;
69  flowRate_ = Function1<scalar>::New("volumetricFlowRate", dict);
70  rhoName_ = "rho";
71  }
72  else if (dict.found("massFlowRate"))
73  {
74  volumetric_ = false;
75  flowRate_ = Function1<scalar>::New("massFlowRate", dict);
76  rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
77  }
78  else
79  {
81  (
82  dict
83  ) << "Please supply either 'volumetricFlowRate' or"
84  << " 'massFlowRate' and 'rho'" << exit(FatalIOError);
85  }
86 
87  // Value field require if mass based
88  if (dict.found("value"))
89  {
91  (
92  vectorField("value", dict, p.size())
93  );
94  }
95  else
96  {
97  evaluate(Pstream::commsTypes::blocking);
98  }
99 }
100 
101 
104 (
106  const fvPatch& p,
108  const fvPatchFieldMapper& mapper
109 )
110 :
111  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
112  flowRate_(ptf.flowRate_, false),
113  volumetric_(ptf.volumetric_),
114  rhoName_(ptf.rhoName_),
115  rhoInlet_(ptf.rhoInlet_),
116  alphaName_(ptf.alphaName_),
117  extrapolateProfile_(ptf.extrapolateProfile_)
118 {}
119 
120 
123 (
125 )
126 :
128  flowRate_(ptf.flowRate_, false),
129  volumetric_(ptf.volumetric_),
130  rhoName_(ptf.rhoName_),
131  rhoInlet_(ptf.rhoInlet_),
132  alphaName_(ptf.alphaName_),
133  extrapolateProfile_(ptf.extrapolateProfile_)
134 {}
135 
136 
139 (
142 )
143 :
145  flowRate_(ptf.flowRate_, false),
146  volumetric_(ptf.volumetric_),
147  rhoName_(ptf.rhoName_),
148  rhoInlet_(ptf.rhoInlet_),
149  alphaName_(ptf.alphaName_),
150  extrapolateProfile_(ptf.extrapolateProfile_)
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
156 template<class AlphaType, class RhoType>
157 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
158 (
159  const AlphaType& alpha,
160  const RhoType& rho
161 )
162 {
163  const scalar t = db().time().timeOutputValue();
164 
165  const vectorField n(patch().nf());
166 
167  if (extrapolateProfile_)
168  {
169  vectorField Up(this->patchInternalField());
170 
171  // Patch normal extrapolated velocity
172  scalarField nUp(n & Up);
173 
174  // Remove the normal component of the extrapolate patch velocity
175  Up -= nUp*n;
176 
177  // Remove any reverse flow
178  nUp = min(nUp, scalar(0));
179 
180  const scalar flowRate = flowRate_->value(t);
181  const scalar estimatedFlowRate =
182  -gSum(alpha*rho*(this->patch().magSf()*nUp));
183 
184  if (estimatedFlowRate/flowRate > 0.5)
185  {
186  nUp *= mag(flowRate)/mag(estimatedFlowRate);
187  }
188  else
189  {
190  nUp -=
191  (flowRate - estimatedFlowRate)
192  /gSum(alpha*rho*patch().magSf());
193  }
194 
195  // Add the corrected normal component of velocity to the patch velocity
196  Up += nUp*n;
197 
198  // Correct the patch velocity
199  this->operator==(Up);
200  }
201  else
202  {
203  const scalar avgU =
204  -flowRate_->value(t)/gSum(alpha*rho*patch().magSf());
205  operator==(avgU*n);
206  }
207 }
208 
209 
210 template<class AlphaType>
211 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
212 (
213  const AlphaType& alpha
214 )
215 {
216  if (volumetric_ || rhoName_ == "none")
217  {
218  updateValues(alpha, one());
219  }
220  else
221  {
222  // Mass flow-rate
223  if (db().foundObject<volScalarField>(rhoName_))
224  {
225  const fvPatchField<scalar>& rhop =
226  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
227 
228  updateValues(alpha, rhop);
229  }
230  else
231  {
232  // Use constant density
233  if (rhoInlet_ < 0)
234  {
236  << "Did not find registered density field " << rhoName_
237  << " and no constant density 'rhoInlet' specified"
238  << exit(FatalError);
239  }
240 
241  updateValues(alpha, rhoInlet_);
242  }
243  }
244 }
245 
246 
248 {
249  if (updated())
250  {
251  return;
252  }
253 
254  if (alphaName_ != word::null)
255  {
256  const fvPatchField<scalar>& alphap =
257  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
258 
259  updateValues(alphap);
260  }
261  else
262  {
263  updateValues(one());
264  }
265 
266  fixedValueFvPatchVectorField::updateCoeffs();
267 }
268 
269 
271 {
273  writeEntry(os, flowRate_());
274  if (!volumetric_)
275  {
276  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
277  writeEntryIfDifferent<scalar>(os, "rhoInlet", -vGreat, rhoInlet_);
278  }
279  writeEntryIfDifferent<word>(os, "alpha", word::null, alphaName_);
280  writeEntry(os, "extrapolateProfile", extrapolateProfile_);
281  writeEntry(os, "value", *this);
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 namespace Foam
288 {
290  (
293  );
294 }
295 
296 
297 // ************************************************************************* //
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Type gSum(const FieldField< Field, Type > &f)
virtual Type value(const scalar x) const =0
Return value as a function of (scalar) independent variable.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
flowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static const word null
An empty word.
Definition: word.H:77
Velocity inlet boundary condition either correcting the extrapolated velocity or creating a uniform v...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:50
IOerror FatalIOError