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-2023 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  const dictionary& dict
39 )
40 :
41  fixedValueFvPatchField<vector>(p, iF, dict, false),
42  rhoOutlet_(dict.lookupOrDefault<scalar>("rhoOutlet", -vGreat))
43 {
44  if (dict.found("volumetricFlowRate"))
45  {
46  volumetric_ = true;
47  flowRate_ = Function1<scalar>::New("volumetricFlowRate", dict);
48  rhoName_ = "rho";
49  }
50  else if (dict.found("massFlowRate"))
51  {
52  volumetric_ = false;
53  flowRate_ = Function1<scalar>::New("massFlowRate", dict);
54  rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
55  }
56  else
57  {
59  (
60  dict
61  ) << "Please supply either 'volumetricFlowRate' or"
62  << " 'massFlowRate' and 'rho'" << exit(FatalIOError);
63  }
64 
65  // Value field require if mass based
66  if (dict.found("value"))
67  {
69  (
70  vectorField("value", dict, p.size())
71  );
72  }
73  else
74  {
76  }
77 }
78 
79 
82 (
84  const fvPatch& p,
86  const fvPatchFieldMapper& mapper
87 )
88 :
89  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
90  flowRate_(ptf.flowRate_, false),
91  volumetric_(ptf.volumetric_),
92  rhoName_(ptf.rhoName_),
93  rhoOutlet_(ptf.rhoOutlet_)
94 {}
95 
96 
99 (
102 )
103 :
104  fixedValueFvPatchField<vector>(ptf, iF),
105  flowRate_(ptf.flowRate_, false),
106  volumetric_(ptf.volumetric_),
107  rhoName_(ptf.rhoName_),
108  rhoOutlet_(ptf.rhoOutlet_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 template<class RhoType>
115 void Foam::flowRateOutletVelocityFvPatchVectorField::updateValues
116 (
117  const RhoType& rho
118 )
119 {
120  const scalar t = db().time().userTimeValue();
121 
122  const vectorField n(patch().nf());
123 
124  // Extrapolate patch velocity
125  vectorField Up(this->patchInternalField());
126 
127  // Patch normal extrapolated velocity
128  scalarField nUp(n & Up);
129 
130  // Remove the normal component of the extrapolate patch velocity
131  Up -= nUp*n;
132 
133  // Remove any reverse flow
134  nUp = max(nUp, scalar(0));
135 
136  const scalar flowRate = flowRate_->value(t);
137  const scalar estimatedFlowRate = gSum(rho*(this->patch().magSf()*nUp));
138 
139  if (estimatedFlowRate/flowRate > 0.5)
140  {
141  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
142  }
143  else
144  {
145  nUp += ((flowRate - estimatedFlowRate)/gSum(rho*patch().magSf()));
146  }
147 
148  // Add the corrected normal component of velocity to the patch velocity
149  Up += nUp*n;
150 
151  // Correct the patch velocity
152  this->operator==(Up);
153 }
154 
155 
157 {
158  if (updated())
159  {
160  return;
161  }
162 
163  if (volumetric_ || rhoName_ == "none")
164  {
165  updateValues(one());
166  }
167  else
168  {
169  // Mass flow-rate
170  if (db().foundObject<volScalarField>(rhoName_))
171  {
172  const fvPatchField<scalar>& rhop =
173  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
174 
175  updateValues(rhop);
176  }
177  else
178  {
179  // Use constant density
180  if (rhoOutlet_ < 0)
181  {
183  << "Did not find registered density field " << rhoName_
184  << " and no constant density 'rhoOutlet' specified"
185  << exit(FatalError);
186  }
187 
188  updateValues(rhoOutlet_);
189  }
190  }
191 
192  fixedValueFvPatchVectorField::updateCoeffs();
193 }
194 
195 
197 {
199  writeEntry(os, flowRate_());
200  if (!volumetric_)
201  {
202  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
203  writeEntryIfDifferent<scalar>(os, "rhoOutlet", -vGreat, rhoOutlet_);
204  }
205  writeEntry(os, "value", *this);
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 namespace Foam
212 {
214  (
217  );
218 }
219 
220 
221 // ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
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.
flowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:51
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p