matchedFlowRateOutletVelocityFvPatchVectorField.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  inletPatchName_(dict.lookup("inletPatch")),
43  volumetric_(dict.lookupOrDefault("volumetric", true))
44 {
45  if (volumetric_)
46  {
47  rhoName_ = "none";
48  }
49  else
50  {
51  rhoName_ = word(dict.lookupOrDefault<word>("rho", "rho"));
52  }
53 
54  // Value field require if mass based
55  if (dict.found("value"))
56  {
58  (
59  vectorField("value", dict, p.size())
60  );
61  }
62  else
63  {
65  }
66 }
67 
68 
71 (
73  const fvPatch& p,
75  const fvPatchFieldMapper& mapper
76 )
77 :
78  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
79  inletPatchName_(ptf.inletPatchName_),
80  volumetric_(ptf.volumetric_),
81  rhoName_(ptf.rhoName_)
82 {}
83 
84 
87 (
90 )
91 :
93  inletPatchName_(ptf.inletPatchName_),
94  volumetric_(ptf.volumetric_),
95  rhoName_(ptf.rhoName_)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
101 template<class RhoType>
102 void Foam::matchedFlowRateOutletVelocityFvPatchVectorField::updateValues
103 (
104  const label inletPatchID,
105  const RhoType& rhoOutlet,
106  const RhoType& rhoInlet
107 )
108 {
109  const fvPatch& p = patch();
110  const fvPatch& inletPatch = p.boundaryMesh()[inletPatchID];
111 
112  const vectorField n(p.nf());
113 
114  // Extrapolate patch velocity
115  vectorField Up(patchInternalField());
116 
117  // Patch normal extrapolated velocity
118  scalarField nUp(n & Up);
119 
120  // Remove the normal component of the extrapolate patch velocity
121  Up -= nUp*n;
122 
123  // Remove any reverse flow
124  nUp = max(nUp, scalar(0));
125 
126  // Lookup non-const access to velocity field
128  (
129  const_cast<volVectorField&>
130  (
131  dynamic_cast<const volVectorField&>(internalField())
132  )
133  );
134 
135  // Get the corresponding inlet velocity patch field
136  fvPatchVectorField& inletPatchU = U.boundaryFieldRef()[inletPatchID];
137 
138  // Ensure that the corresponding inlet velocity patch field is up-to-date
139  inletPatchU.updateCoeffs();
140 
141  // Calculate the inlet patch flow rate
142  const scalar flowRate = -gSum(rhoInlet*(inletPatch.Sf() & inletPatchU));
143 
144  // Calculate the extrapolated outlet patch flow rate
145  const scalar estimatedFlowRate = gSum(rhoOutlet*(patch().magSf()*nUp));
146 
147  if (estimatedFlowRate/flowRate > 0.5)
148  {
149  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
150  }
151  else
152  {
153  nUp += ((flowRate - estimatedFlowRate)/gSum(rhoOutlet*patch().magSf()));
154  }
155 
156  // Add the corrected normal component of velocity to the patch velocity
157  Up += nUp*n;
158 
159  // Correct the patch velocity
160  operator==(Up);
161 }
162 
163 
165 {
166  if (updated())
167  {
168  return;
169  }
170 
171  // Find corresponding inlet patch
172  const label inletPatchID =
173  patch().patch().boundaryMesh().findPatchID(inletPatchName_);
174 
175  if (inletPatchID < 0)
176  {
178  << "Unable to find inlet patch " << inletPatchName_
179  << exit(FatalError);
180  }
181 
182  if (volumetric_)
183  {
184  updateValues(inletPatchID, one(), one());
185  }
186  else
187  {
188  // Mass flow-rate
189  if (db().foundObject<volScalarField>(rhoName_))
190  {
191  const volScalarField& rho = db().lookupObject<volScalarField>
192  (
193  rhoName_
194  );
195 
196  updateValues
197  (
198  inletPatchID,
199  rho.boundaryField()[patch().index()],
200  rho.boundaryField()[inletPatchID]
201  );
202  }
203  else
204  {
206  << "Cannot find density field " << rhoName_ << exit(FatalError);
207  }
208  }
209 
210  fixedValueFvPatchVectorField::updateCoeffs();
211 }
212 
213 
215 (
216  Ostream& os
217 ) const
218 {
220  writeEntry(os, "inletPatch", inletPatchName_);
221  if (!volumetric_)
222  {
223  writeEntry(os, "volumetric", volumetric_);
224  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
225  }
226  writeEntry(os, "value", *this);
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
231 
232 namespace Foam
233 {
235  (
238  );
239 }
240 
241 
242 // ************************************************************************* //
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...
Generic GeometricField class.
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...
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:136
Velocity outlet boundary condition which corrects the extrapolated velocity to match the flow rate of...
matchedFlowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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 FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Type gSum(const FieldField< Field, Type > &f)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
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