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 (
126 )
127 :
129  flowRate_(ptf.flowRate_, false),
130  volumetric_(ptf.volumetric_),
131  rhoName_(ptf.rhoName_),
132  rhoInlet_(ptf.rhoInlet_),
133  alphaName_(ptf.alphaName_),
134  extrapolateProfile_(ptf.extrapolateProfile_)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
140 template<class AlphaType, class RhoType>
141 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
142 (
143  const AlphaType& alpha,
144  const RhoType& rho
145 )
146 {
147  const scalar t = db().time().timeOutputValue();
148 
149  const vectorField n(patch().nf());
150 
151  if (extrapolateProfile_)
152  {
153  vectorField Up(this->patchInternalField());
154 
155  // Patch normal extrapolated velocity
156  scalarField nUp(n & Up);
157 
158  // Remove the normal component of the extrapolate patch velocity
159  Up -= nUp*n;
160 
161  // Remove any reverse flow
162  nUp = min(nUp, scalar(0));
163 
164  const scalar flowRate = flowRate_->value(t);
165  const scalar estimatedFlowRate =
166  -gSum(alpha*rho*(this->patch().magSf()*nUp));
167 
168  if (estimatedFlowRate/flowRate > 0.5)
169  {
170  nUp *= mag(flowRate)/mag(estimatedFlowRate);
171  }
172  else
173  {
174  nUp -=
175  (flowRate - estimatedFlowRate)
176  /gSum(alpha*rho*patch().magSf());
177  }
178 
179  // Add the corrected normal component of velocity to the patch velocity
180  Up += nUp*n;
181 
182  // Correct the patch velocity
183  this->operator==(Up);
184  }
185  else
186  {
187  const scalar avgU =
188  -flowRate_->value(t)/gSum(alpha*rho*patch().magSf());
189  operator==(avgU*n);
190  }
191 }
192 
193 
194 template<class AlphaType>
195 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
196 (
197  const AlphaType& alpha
198 )
199 {
200  if (volumetric_ || rhoName_ == "none")
201  {
202  updateValues(alpha, one());
203  }
204  else
205  {
206  // Mass flow-rate
207  if (db().foundObject<volScalarField>(rhoName_))
208  {
209  const fvPatchField<scalar>& rhop =
210  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
211 
212  updateValues(alpha, rhop);
213  }
214  else
215  {
216  // Use constant density
217  if (rhoInlet_ < 0)
218  {
220  << "Did not find registered density field " << rhoName_
221  << " and no constant density 'rhoInlet' specified"
222  << exit(FatalError);
223  }
224 
225  updateValues(alpha, rhoInlet_);
226  }
227  }
228 }
229 
230 
232 {
233  if (updated())
234  {
235  return;
236  }
237 
238  if (alphaName_ != word::null)
239  {
240  const fvPatchField<scalar>& alphap =
241  patch().lookupPatchField<volScalarField, scalar>(alphaName_);
242 
243  updateValues(alphap);
244  }
245  else
246  {
247  updateValues(one());
248  }
249 
250  fixedValueFvPatchVectorField::updateCoeffs();
251 }
252 
253 
255 {
257  writeEntry(os, flowRate_());
258  if (!volumetric_)
259  {
260  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
261  writeEntryIfDifferent<scalar>(os, "rhoInlet", -vGreat, rhoInlet_);
262  }
263  writeEntryIfDifferent<word>(os, "alpha", word::null, alphaName_);
264  writeEntry(os, "extrapolateProfile", extrapolateProfile_);
265  writeEntry(os, "value", *this);
266 }
267 
268 
269 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
270 
271 namespace Foam
272 {
274  (
277  );
278 }
279 
280 
281 // ************************************************************************* //
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
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
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:62
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:260
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 x.
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:335
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