flowRateInletVelocityFvPatchVectorField.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-2017 OpenFOAM Foundation
9  Copyright (C) 2019-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 "one.H"
33 #include "Switch.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
45  flowRate_(nullptr),
46  rhoName_("rho"),
47  rhoInlet_(0),
48  volumetric_(false),
49  extrapolateProfile_(false)
50 {}
51 
52 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
62  flowRate_(nullptr),
63  rhoName_("rho"),
64  rhoInlet_(dict.getOrDefault<scalar>("rhoInlet", -VGREAT)),
65  volumetric_(false),
66  extrapolateProfile_
67  (
68  dict.getOrDefault<Switch>("extrapolateProfile", false)
69  )
70 {
71  flowRate_ =
72  Function1<scalar>::NewIfPresent("volumetricFlowRate", dict, &db());
73 
74  if (flowRate_)
75  {
76  volumetric_ = true;
77  }
78  else
79  {
80  dict.readIfPresent("rho", rhoName_);
81  flowRate_ =
82  Function1<scalar>::NewIfPresent("massFlowRate", dict, &db());
83  }
84 
85  if (!flowRate_)
86  {
88  << "Please supply either 'volumetricFlowRate' or"
89  << " 'massFlowRate' (optional: with 'rho')" << nl
90  << exit(FatalIOError);
91  }
92 
93  // Value field required if mass based
94  if (!this->readValueEntry(dict))
95  {
97  }
98 }
99 
100 
103 (
104  const flowRateInletVelocityFvPatchVectorField& ptf,
105  const fvPatch& p,
106  const DimensionedField<vector, volMesh>& iF,
107  const fvPatchFieldMapper& mapper
108 )
109 :
110  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
111  flowRate_(ptf.flowRate_.clone()),
112  rhoName_(ptf.rhoName_),
113  rhoInlet_(ptf.rhoInlet_),
114  volumetric_(ptf.volumetric_),
115  extrapolateProfile_(ptf.extrapolateProfile_)
116 {}
117 
118 
121 (
123 )
124 :
126  flowRate_(ptf.flowRate_.clone()),
127  rhoName_(ptf.rhoName_),
128  rhoInlet_(ptf.rhoInlet_),
129  volumetric_(ptf.volumetric_),
130  extrapolateProfile_(ptf.extrapolateProfile_)
131 {}
132 
133 
136 (
139 )
140 :
141  fixedValueFvPatchField<vector>(ptf, iF),
142  flowRate_(ptf.flowRate_.clone()),
143  rhoName_(ptf.rhoName_),
144  rhoInlet_(ptf.rhoInlet_),
145  volumetric_(ptf.volumetric_),
146  extrapolateProfile_(ptf.extrapolateProfile_)
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
152 template<class RhoType>
153 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
154 (
155  const RhoType& rho
156 )
157 {
158  const scalar t = db().time().timeOutputValue();
159 
160  const vectorField n(patch().nf());
161 
162  if (extrapolateProfile_)
163  {
164  vectorField Up(this->patchInternalField());
165 
166  // Patch normal extrapolated velocity
167  scalarField nUp(n & Up);
168 
169  // Remove the normal component of the extrapolate patch velocity
170  Up -= nUp*n;
171 
172  // Remove any reverse flow
173  nUp = min(nUp, scalar(0));
174 
175  const scalar flowRate = flowRate_->value(t);
176  const scalar estimatedFlowRate = -gSum(rho*(this->patch().magSf()*nUp));
177 
178  if (estimatedFlowRate > 0.5*flowRate)
179  {
180  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
181  }
182  else
183  {
184  nUp -= ((flowRate - estimatedFlowRate)/gSum(rho*patch().magSf()));
185  }
186 
187  // Add the corrected normal component of velocity to the patch velocity
188  Up += nUp*n;
189 
190  // Correct the patch velocity
191  this->operator==(Up);
192  }
193  else
194  {
195  const scalar avgU = -flowRate_->value(t)/gSum(rho*patch().magSf());
196  operator==(avgU*n);
197  }
198 }
199 
200 
202 {
203  if (updated())
204  {
205  return;
206  }
207 
208  if (volumetric_ || rhoName_ == "none")
209  {
210  updateValues(one{});
211  }
212  else
213  {
214  // Mass flow-rate
215  if (db().foundObject<volScalarField>(rhoName_))
216  {
217  const auto& rhop =
218  patch().lookupPatchField<volScalarField>(rhoName_);
219 
220  updateValues(rhop);
221  }
222  else
223  {
224  // Use constant density
225  if (rhoInlet_ < 0)
226  {
228  << "Did not find registered density field " << rhoName_
229  << " and no constant density 'rhoInlet' specified"
230  << exit(FatalError);
231  }
232 
233  updateValues(rhoInlet_);
234  }
235  }
236 
237  fixedValueFvPatchVectorField::updateCoeffs();
238 }
239 
240 
242 {
244  flowRate_->writeData(os);
245  if (!volumetric_)
246  {
247  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
248  os.writeEntryIfDifferent<scalar>("rhoInlet", -VGREAT, rhoInlet_);
249  }
250  if (extrapolateProfile_)
251  {
252  os.writeEntry("extrapolateProfile", extrapolateProfile_);
253  }
255 }
256 
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 namespace Foam
261 {
263  (
265  flowRateInletVelocityFvPatchVectorField
266  );
267 }
268 
269 
270 // ************************************************************************* //
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
"blocking" : (MPI_Bsend, MPI_Recv)
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Type gSum(const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
flowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Velocity inlet boundary condition either correcting the extrapolated velocity or creating a uniform v...
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
static autoPtr< Function1< Type > > NewIfPresent(const word &entryName, const dictionary &dict, const word &redirectType, const objectRegistry *obrPtr=nullptr)
An optional selector, with fallback redirection.
Definition: Function1New.C:209
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...