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 | Copyright (C) 2011-2017 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  extrapolateProfile_(false)
46 {}
47 
48 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
58  rhoInlet_(dict.lookupOrDefault<scalar>("rhoInlet", -VGREAT)),
59  extrapolateProfile_
60  (
61  dict.lookupOrDefault<Switch>("extrapolateProfile", false)
62  )
63 {
64  if (dict.found("volumetricFlowRate"))
65  {
66  volumetric_ = true;
67  flowRate_ = Function1<scalar>::New("volumetricFlowRate", dict);
68  rhoName_ = "rho";
69  }
70  else if (dict.found("massFlowRate"))
71  {
72  volumetric_ = false;
73  flowRate_ = Function1<scalar>::New("massFlowRate", dict);
74  rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
75  }
76  else
77  {
79  (
80  dict
81  ) << "Please supply either 'volumetricFlowRate' or"
82  << " 'massFlowRate' and 'rho'" << exit(FatalIOError);
83  }
84 
85  // Value field require if mass based
86  if (dict.found("value"))
87  {
89  (
90  vectorField("value", dict, p.size())
91  );
92  }
93  else
94  {
95  evaluate(Pstream::commsTypes::blocking);
96  }
97 }
98 
99 
102 (
104  const fvPatch& p,
106  const fvPatchFieldMapper& mapper
107 )
108 :
109  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
110  flowRate_(ptf.flowRate_, false),
111  volumetric_(ptf.volumetric_),
112  rhoName_(ptf.rhoName_),
113  rhoInlet_(ptf.rhoInlet_),
114  extrapolateProfile_(ptf.extrapolateProfile_)
115 {}
116 
117 
120 (
122 )
123 :
125  flowRate_(ptf.flowRate_, false),
126  volumetric_(ptf.volumetric_),
127  rhoName_(ptf.rhoName_),
128  rhoInlet_(ptf.rhoInlet_),
129  extrapolateProfile_(ptf.extrapolateProfile_)
130 {}
131 
132 
135 (
138 )
139 :
141  flowRate_(ptf.flowRate_, false),
142  volumetric_(ptf.volumetric_),
143  rhoName_(ptf.rhoName_),
144  rhoInlet_(ptf.rhoInlet_),
145  extrapolateProfile_(ptf.extrapolateProfile_)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 template<class RhoType>
152 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
153 (
154  const RhoType& rho
155 )
156 {
157  const scalar t = db().time().timeOutputValue();
158 
159  const vectorField n(patch().nf());
160 
161  if (extrapolateProfile_)
162  {
163  vectorField Up(this->patchInternalField());
164 
165  // Patch normal extrapolated velocity
166  scalarField nUp(n & Up);
167 
168  // Remove the normal component of the extrapolate patch velocity
169  Up -= nUp*n;
170 
171  // Remove any reverse flow
172  nUp = min(nUp, scalar(0));
173 
174  const scalar flowRate = flowRate_->value(t);
175  const scalar estimatedFlowRate = -gSum(rho*(this->patch().magSf()*nUp));
176 
177  if (estimatedFlowRate/flowRate > 0.5)
178  {
179  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
180  }
181  else
182  {
183  nUp -= ((flowRate - estimatedFlowRate)/gSum(rho*patch().magSf()));
184  }
185 
186  // Add the corrected normal component of velocity to the patch velocity
187  Up += nUp*n;
188 
189  // Correct the patch velocity
190  this->operator==(Up);
191  }
192  else
193  {
194  const scalar avgU = -flowRate_->value(t)/gSum(rho*patch().magSf());
195  operator==(avgU*n);
196  }
197 }
198 
199 
201 {
202  if (updated())
203  {
204  return;
205  }
206 
207  if (volumetric_ || rhoName_ == "none")
208  {
209  updateValues(one());
210  }
211  else
212  {
213  // Mass flow-rate
214  if (db().foundObject<volScalarField>(rhoName_))
215  {
216  const fvPatchField<scalar>& rhop =
217  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
218 
219  updateValues(rhop);
220  }
221  else
222  {
223  // Use constant density
224  if (rhoInlet_ < 0)
225  {
227  << "Did not find registered density field " << rhoName_
228  << " and no constant density 'rhoInlet' specified"
229  << exit(FatalError);
230  }
231 
232  updateValues(rhoInlet_);
233  }
234  }
235 
236  fixedValueFvPatchVectorField::updateCoeffs();
237 }
238 
239 
241 {
243  flowRate_->writeData(os);
244  if (!volumetric_)
245  {
246  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
247  writeEntryIfDifferent<scalar>(os, "rhoInlet", -VGREAT, rhoInlet_);
248  }
249  os.writeKeyword("extrapolateProfile")
250  << extrapolateProfile_ << token::END_STATEMENT << nl;
251  writeEntry("value", os);
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 namespace Foam
258 {
260  (
263  );
264 }
265 
266 
267 // ************************************************************************* //
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:137
#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
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.
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:362
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
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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
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.
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:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
#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
volScalarField & p
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