flowRateOutletVelocityFvPatchVectorField.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) 2017-2019 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 
27 #include "volFields.H"
28 #include "one.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38 )
39 :
41  flowRate_(),
42  volumetric_(false),
43  rhoName_("rho"),
44  rhoOutlet_(0.0)
45 {}
46 
47 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
57  rhoOutlet_(dict.lookupOrDefault<scalar>("rhoOutlet", -vGreat))
58 {
59  if (dict.found("volumetricFlowRate"))
60  {
61  volumetric_ = true;
62  flowRate_ = Function1<scalar>::New("volumetricFlowRate", dict);
63  rhoName_ = "rho";
64  }
65  else if (dict.found("massFlowRate"))
66  {
67  volumetric_ = false;
68  flowRate_ = Function1<scalar>::New("massFlowRate", dict);
69  rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
70  }
71  else
72  {
74  (
75  dict
76  ) << "Please supply either 'volumetricFlowRate' or"
77  << " 'massFlowRate' and 'rho'" << exit(FatalIOError);
78  }
79 
80  // Value field require if mass based
81  if (dict.found("value"))
82  {
84  (
85  vectorField("value", dict, p.size())
86  );
87  }
88  else
89  {
90  evaluate(Pstream::commsTypes::blocking);
91  }
92 }
93 
94 
97 (
99  const fvPatch& p,
101  const fvPatchFieldMapper& mapper
102 )
103 :
104  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
105  flowRate_(ptf.flowRate_, false),
106  volumetric_(ptf.volumetric_),
107  rhoName_(ptf.rhoName_),
108  rhoOutlet_(ptf.rhoOutlet_)
109 {}
110 
111 
114 (
116 )
117 :
119  flowRate_(ptf.flowRate_, false),
120  volumetric_(ptf.volumetric_),
121  rhoName_(ptf.rhoName_),
122  rhoOutlet_(ptf.rhoOutlet_)
123 {}
124 
125 
128 (
131 )
132 :
134  flowRate_(ptf.flowRate_, false),
135  volumetric_(ptf.volumetric_),
136  rhoName_(ptf.rhoName_),
137  rhoOutlet_(ptf.rhoOutlet_)
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/flowRate > 0.5)
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 fvPatchField<scalar>& rhop =
202  patch().lookupPatchField<volScalarField, scalar>(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  writeEntry(os, flowRate_());
229  if (!volumetric_)
230  {
231  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
232  writeEntryIfDifferent<scalar>(os, "rhoOutlet", -vGreat, rhoOutlet_);
233  }
234  writeEntry(os, "value", *this);
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
243  (
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
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#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:163
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
#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
Velocity outlet boundary condition which corrects the extrapolated velocity to match the specified fl...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
flowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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