flowRateOutletVelocityFvPatchVectorField.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) 2017 OpenFOAM Foundation
9  Copyright (C) 2020-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 
30 #include "volFields.H"
31 #include "one.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
44  flowRate_(nullptr),
45  rhoName_("rho"),
46  rhoOutlet_(0),
47  volumetric_(false)
48 {}
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
60  flowRate_(nullptr),
61  rhoName_("rho"),
62  rhoOutlet_(dict.getOrDefault<scalar>("rhoOutlet", -VGREAT)),
63  volumetric_(false)
64 {
65  flowRate_ =
66  Function1<scalar>::NewIfPresent("volumetricFlowRate", dict, &db());
67 
68  if (flowRate_)
69  {
70  volumetric_ = true;
71  }
72  else
73  {
74  dict.readIfPresent("rho", rhoName_);
75  flowRate_ =
76  Function1<scalar>::NewIfPresent("massFlowRate", dict, &db());
77  }
78 
79  if (!flowRate_)
80  {
82  << "Please supply either 'volumetricFlowRate' or"
83  << " 'massFlowRate' (optional: with 'rho')" << nl
84  << exit(FatalIOError);
85  }
86 
87  // Value field required if mass based
88  if (!this->readValueEntry(dict))
89  {
91  }
92 }
93 
94 
97 (
98  const flowRateOutletVelocityFvPatchVectorField& ptf,
99  const fvPatch& p,
100  const DimensionedField<vector, volMesh>& iF,
101  const fvPatchFieldMapper& mapper
102 )
103 :
104  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
105  flowRate_(ptf.flowRate_.clone()),
106  rhoName_(ptf.rhoName_),
107  rhoOutlet_(ptf.rhoOutlet_),
108  volumetric_(ptf.volumetric_)
109 {}
110 
111 
114 (
116 )
117 :
119  flowRate_(ptf.flowRate_.clone()),
120  rhoName_(ptf.rhoName_),
121  rhoOutlet_(ptf.rhoOutlet_),
122  volumetric_(ptf.volumetric_)
123 {}
124 
125 
128 (
131 )
132 :
133  fixedValueFvPatchField<vector>(ptf, iF),
134  flowRate_(ptf.flowRate_.clone()),
135  rhoName_(ptf.rhoName_),
136  rhoOutlet_(ptf.rhoOutlet_),
137  volumetric_(ptf.volumetric_)
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
143 template<class RhoType>
144 void Foam::flowRateOutletVelocityFvPatchVectorField::updateValues
145 (
146  const RhoType& rho
147 )
148 {
149  const scalar t = db().time().timeOutputValue();
150 
151  const vectorField n(patch().nf());
152 
153  // Extrapolate patch velocity
154  vectorField Up(this->patchInternalField());
155 
156  // Patch normal extrapolated velocity
157  scalarField nUp(n & Up);
158 
159  // Remove the normal component of the extrapolate patch velocity
160  Up -= nUp*n;
161 
162  // Remove any reverse flow
163  nUp = max(nUp, scalar(0));
164 
165  const scalar flowRate = flowRate_->value(t);
166  const scalar estimatedFlowRate = gSum(rho*(this->patch().magSf()*nUp));
167 
168  if (estimatedFlowRate > 0.5*flowRate)
169  {
170  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
171  }
172  else
173  {
174  nUp += ((flowRate - estimatedFlowRate)/gSum(rho*patch().magSf()));
175  }
176 
177  // Add the corrected normal component of velocity to the patch velocity
178  Up += nUp*n;
179 
180  // Correct the patch velocity
181  this->operator==(Up);
182 }
183 
184 
186 {
187  if (updated())
188  {
189  return;
190  }
191 
192  if (volumetric_ || rhoName_ == "none")
193  {
194  updateValues(one{});
195  }
196  else
197  {
198  // Mass flow-rate
199  if (db().foundObject<volScalarField>(rhoName_))
200  {
201  const auto& rhop =
202  patch().lookupPatchField<volScalarField>(rhoName_);
203 
204  updateValues(rhop);
205  }
206  else
207  {
208  // Use constant density
209  if (rhoOutlet_ < 0)
210  {
212  << "Did not find registered density field " << rhoName_
213  << " and no constant density 'rhoOutlet' specified"
214  << exit(FatalError);
215  }
216 
217  updateValues(rhoOutlet_);
218  }
219  }
220 
221  fixedValueFvPatchVectorField::updateCoeffs();
222 }
223 
224 
226 {
228  flowRate_->writeData(os);
229  if (!volumetric_)
230  {
231  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
232  os.writeEntryIfDifferent<scalar>("rhoOutlet", -VGREAT, rhoOutlet_);
233  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
243  (
245  flowRateOutletVelocityFvPatchVectorField
246  );
247 }
248 
249 
250 // ************************************************************************* //
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
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
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
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.
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.
Vector< scalar > vector
Definition: vector.H:57
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
Velocity outlet boundary condition which corrects the extrapolated velocity to match the specified fl...
Field< vector > vectorField
Specialisation of Field<T> for vector.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
flowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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 ...