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-2020 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 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
44  flowRate_(),
45  rhoName_("rho"),
46  rhoInlet_(0.0),
47  volumetric_(false),
48  extrapolateProfile_(false)
49 {}
50 
51 
54 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
61  rhoName_("rho"),
62  rhoInlet_(dict.getOrDefault<scalar>("rhoInlet", -VGREAT)),
63  volumetric_(false),
64  extrapolateProfile_
65  (
66  dict.getOrDefault<Switch>("extrapolateProfile", false)
67  )
68 {
69  if (dict.found("volumetricFlowRate"))
70  {
71  volumetric_ = true;
72  flowRate_ = Function1<scalar>::New("volumetricFlowRate", dict);
73  }
74  else if (dict.found("massFlowRate"))
75  {
76  volumetric_ = false;
77  flowRate_ = Function1<scalar>::New("massFlowRate", dict);
78  rhoName_ = dict.getOrDefault<word>("rho", "rho");
79  }
80  else
81  {
83  << "Please supply either 'volumetricFlowRate' or"
84  << " 'massFlowRate' and 'rho'" << nl
85  << exit(FatalIOError);
86  }
87 
88  // Value field require if mass based
89  if (dict.found("value"))
90  {
92  (
93  vectorField("value", dict, p.size())
94  );
95  }
96  else
97  {
98  evaluate(Pstream::commsTypes::blocking);
99  }
100 }
101 
102 
105 (
107  const fvPatch& p,
109  const fvPatchFieldMapper& mapper
110 )
111 :
112  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
113  flowRate_(ptf.flowRate_.clone()),
114  rhoName_(ptf.rhoName_),
115  rhoInlet_(ptf.rhoInlet_),
116  volumetric_(ptf.volumetric_),
117  extrapolateProfile_(ptf.extrapolateProfile_)
118 {}
119 
120 
123 (
125 )
126 :
128  flowRate_(ptf.flowRate_.clone()),
129  rhoName_(ptf.rhoName_),
130  rhoInlet_(ptf.rhoInlet_),
131  volumetric_(ptf.volumetric_),
132  extrapolateProfile_(ptf.extrapolateProfile_)
133 {}
134 
135 
138 (
141 )
142 :
144  flowRate_(ptf.flowRate_.clone()),
145  rhoName_(ptf.rhoName_),
146  rhoInlet_(ptf.rhoInlet_),
147  volumetric_(ptf.volumetric_),
148  extrapolateProfile_(ptf.extrapolateProfile_)
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 template<class RhoType>
155 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
156 (
157  const RhoType& rho
158 )
159 {
160  const scalar t = db().time().timeOutputValue();
161 
162  const vectorField n(patch().nf());
163 
164  if (extrapolateProfile_)
165  {
166  vectorField Up(this->patchInternalField());
167 
168  // Patch normal extrapolated velocity
169  scalarField nUp(n & Up);
170 
171  // Remove the normal component of the extrapolate patch velocity
172  Up -= nUp*n;
173 
174  // Remove any reverse flow
175  nUp = min(nUp, scalar(0));
176 
177  const scalar flowRate = flowRate_->value(t);
178  const scalar estimatedFlowRate = -gSum(rho*(this->patch().magSf()*nUp));
179 
180  if (estimatedFlowRate/flowRate > 0.5)
181  {
182  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
183  }
184  else
185  {
186  nUp -= ((flowRate - estimatedFlowRate)/gSum(rho*patch().magSf()));
187  }
188 
189  // Add the corrected normal component of velocity to the patch velocity
190  Up += nUp*n;
191 
192  // Correct the patch velocity
193  this->operator==(Up);
194  }
195  else
196  {
197  const scalar avgU = -flowRate_->value(t)/gSum(rho*patch().magSf());
198  operator==(avgU*n);
199  }
200 }
201 
202 
204 {
205  if (updated())
206  {
207  return;
208  }
209 
210  if (volumetric_ || rhoName_ == "none")
211  {
212  updateValues(one{});
213  }
214  else
215  {
216  // Mass flow-rate
217  if (db().foundObject<volScalarField>(rhoName_))
218  {
219  const fvPatchField<scalar>& rhop =
220  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
221 
222  updateValues(rhop);
223  }
224  else
225  {
226  // Use constant density
227  if (rhoInlet_ < 0)
228  {
230  << "Did not find registered density field " << rhoName_
231  << " and no constant density 'rhoInlet' specified"
232  << exit(FatalError);
233  }
234 
235  updateValues(rhoInlet_);
236  }
237  }
238 
239  fixedValueFvPatchVectorField::updateCoeffs();
240 }
241 
242 
244 {
246  flowRate_->writeData(os);
247  if (!volumetric_)
248  {
249  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
250  os.writeEntryIfDifferent<scalar>("rhoInlet", -VGREAT, rhoInlet_);
251  }
252  os.writeEntry("extrapolateProfile", extrapolateProfile_);
253  writeEntry("value", os);
254 }
255 
256 
257 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258 
259 namespace Foam
260 {
262  (
265  );
266 }
267 
268 
269 // ************************************************************************* //
flowRateInletVelocityFvPatchVectorField.H
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::flowRateInletVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: flowRateInletVelocityFvPatchVectorField.C:243
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
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:244
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
Foam::FatalIOError
IOerror FatalIOError
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:80
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
one.H
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::flowRateInletVelocityFvPatchVectorField::flowRateInletVelocityFvPatchVectorField
flowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: flowRateInletVelocityFvPatchVectorField.C:38
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::flowRateInletVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: flowRateInletVelocityFvPatchVectorField.C:203
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::flowRateInletVelocityFvPatchVectorField
Velocity inlet boundary condition either correcting the extrapolated velocity or creating a uniform v...
Definition: flowRateInletVelocityFvPatchVectorField.H:145
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:401
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::stringOps::evaluate
string evaluate(const std::string &s, size_t pos=0, size_t len=std::string::npos)
Definition: stringOpsEvaluate.C:37
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54