matchedFlowRateOutletVelocityFvPatchVectorField.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  inletPatchName_(),
45  rhoName_("rho"),
46  volumetric_(false)
47 {}
48 
49 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
59  inletPatchName_(dict.get<word>("inletPatch")),
60  rhoName_(),
61  volumetric_(dict.getOrDefault("volumetric", true))
62 {
63  if (volumetric_)
64  {
65  rhoName_ = "none";
66  }
67  else
68  {
69  rhoName_ = dict.getOrDefault<word>("rho", "rho");
70  }
71 
72  // Value field required if mass based
73  if (!this->readValueEntry(dict))
74  {
76  }
77 }
78 
79 
82 (
83  const matchedFlowRateOutletVelocityFvPatchVectorField& ptf,
84  const fvPatch& p,
85  const DimensionedField<vector, volMesh>& iF,
86  const fvPatchFieldMapper& mapper
87 )
88 :
89  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
90  inletPatchName_(ptf.inletPatchName_),
91  rhoName_(ptf.rhoName_),
92  volumetric_(ptf.volumetric_)
93 {}
94 
95 
98 (
100 )
101 :
103  inletPatchName_(ptf.inletPatchName_),
104  rhoName_(ptf.rhoName_),
105  volumetric_(ptf.volumetric_)
106 {}
107 
108 
111 (
114 )
115 :
116  fixedValueFvPatchField<vector>(ptf, iF),
117  inletPatchName_(ptf.inletPatchName_),
118  rhoName_(ptf.rhoName_),
119  volumetric_(ptf.volumetric_)
120 {}
121 
122 
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124 
125 template<class RhoType>
126 void Foam::matchedFlowRateOutletVelocityFvPatchVectorField::updateValues
127 (
128  const label inletPatchID,
129  const RhoType& rhoOutlet,
130  const RhoType& rhoInlet
131 )
132 {
133  const fvPatch& p = patch();
134  const fvPatch& inletPatch = p.boundaryMesh()[inletPatchID];
135 
136  const vectorField n(p.nf());
137 
138  // Extrapolate patch velocity
139  vectorField Up(patchInternalField());
140 
141  // Patch normal extrapolated velocity
142  scalarField nUp(n & Up);
143 
144  // Remove the normal component of the extrapolate patch velocity
145  Up -= nUp*n;
146 
147  // Remove any reverse flow
148  nUp = max(nUp, scalar(0));
149 
150  // Lookup non-const access to velocity field
152  (
153  const_cast<volVectorField&>
154  (
155  dynamic_cast<const volVectorField&>(internalField())
156  )
157  );
158 
159  // Get the corresponding inlet velocity patch field
160  fvPatchVectorField& inletPatchU = U.boundaryFieldRef()[inletPatchID];
161 
162  // Ensure that the corresponding inlet velocity patch field is up-to-date
163  inletPatchU.updateCoeffs();
164 
165  // Calculate the inlet patch flow rate
166  const scalar flowRate = -gSum(rhoInlet*(inletPatch.Sf() & inletPatchU));
167 
168  // Calculate the extrapolated outlet patch flow rate
169  const scalar estimatedFlowRate = gSum(rhoOutlet*(patch().magSf()*nUp));
170 
171  if (estimatedFlowRate > 0.5*flowRate)
172  {
173  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
174  }
175  else
176  {
177  nUp += ((flowRate - estimatedFlowRate)/gSum(rhoOutlet*patch().magSf()));
178  }
179 
180  // Add the corrected normal component of velocity to the patch velocity
181  Up += nUp*n;
182 
183  // Correct the patch velocity
184  operator==(Up);
185 }
186 
187 
189 {
190  if (updated())
191  {
192  return;
193  }
194 
195  // Find corresponding inlet patch
196  const label inletPatchID =
197  patch().patch().boundaryMesh().findPatchID(inletPatchName_);
198 
199  if (inletPatchID < 0)
200  {
202  << "Unable to find inlet patch " << inletPatchName_
203  << exit(FatalError);
204  }
205 
206  if (volumetric_)
207  {
208  updateValues(inletPatchID, one{}, one{});
209  }
210  else
211  {
212  // Mass flow-rate
213  if (db().foundObject<volScalarField>(rhoName_))
214  {
215  const volScalarField& rho = db().lookupObject<volScalarField>
216  (
217  rhoName_
218  );
219 
220  updateValues
221  (
222  inletPatchID,
223  rho.boundaryField()[patch().index()],
224  rho.boundaryField()[inletPatchID]
225  );
226  }
227  else
228  {
230  << "Cannot find density field " << rhoName_ << exit(FatalError);
231  }
232  }
233 
234  fixedValueFvPatchVectorField::updateCoeffs();
235 }
236 
237 
239 (
240  Ostream& os
241 ) const
242 {
244  os.writeEntry("inletPatch", inletPatchName_);
245  if (!volumetric_)
246  {
247  os.writeEntry("volumetric", volumetric_);
248  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
249  }
251 }
252 
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 namespace Foam
257 {
259  (
261  matchedFlowRateOutletVelocityFvPatchVectorField
262  );
263 }
264 
265 
266 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
matchedFlowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Velocity outlet boundary condition which corrects the extrapolated velocity to match the flow rate of...
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
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
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
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)
U
Definition: pEqn.H:72
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
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.