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  volumetric_(false),
46  rhoName_("rho")
47 {}
48 
49 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
59  inletPatchName_(dict.get<word>("inletPatch")),
60  volumetric_(dict.getOrDefault("volumetric", true))
61 {
62  if (volumetric_)
63  {
64  rhoName_ = "none";
65  }
66  else
67  {
68  rhoName_ = dict.getOrDefault<word>("rho", "rho");
69  }
70 
71  // Value field require if mass based
72  if (dict.found("value"))
73  {
75  (
76  vectorField("value", dict, p.size())
77  );
78  }
79  else
80  {
81  evaluate(Pstream::commsTypes::blocking);
82  }
83 }
84 
85 
88 (
90  const fvPatch& p,
92  const fvPatchFieldMapper& mapper
93 )
94 :
95  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
96  inletPatchName_(ptf.inletPatchName_),
97  volumetric_(ptf.volumetric_),
98  rhoName_(ptf.rhoName_)
99 {}
100 
101 
104 (
106 )
107 :
109  inletPatchName_(ptf.inletPatchName_),
110  volumetric_(ptf.volumetric_),
111  rhoName_(ptf.rhoName_)
112 {}
113 
114 
117 (
120 )
121 :
123  inletPatchName_(ptf.inletPatchName_),
124  volumetric_(ptf.volumetric_),
125  rhoName_(ptf.rhoName_)
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
131 template<class RhoType>
132 void Foam::matchedFlowRateOutletVelocityFvPatchVectorField::updateValues
133 (
134  const label inletPatchID,
135  const RhoType& rhoOutlet,
136  const RhoType& rhoInlet
137 )
138 {
139  const fvPatch& p = patch();
140  const fvPatch& inletPatch = p.boundaryMesh()[inletPatchID];
141 
142  const vectorField n(p.nf());
143 
144  // Extrapolate patch velocity
145  vectorField Up(patchInternalField());
146 
147  // Patch normal extrapolated velocity
148  scalarField nUp(n & Up);
149 
150  // Remove the normal component of the extrapolate patch velocity
151  Up -= nUp*n;
152 
153  // Remove any reverse flow
154  nUp = max(nUp, scalar(0));
155 
156  // Lookup non-const access to velocity field
158  (
159  const_cast<volVectorField&>
160  (
161  dynamic_cast<const volVectorField&>(internalField())
162  )
163  );
164 
165  // Get the corresponding inlet velocity patch field
166  fvPatchVectorField& inletPatchU = U.boundaryFieldRef()[inletPatchID];
167 
168  // Ensure that the corresponding inlet velocity patch field is up-to-date
169  inletPatchU.updateCoeffs();
170 
171  // Calculate the inlet patch flow rate
172  const scalar flowRate = -gSum(rhoInlet*(inletPatch.Sf() & inletPatchU));
173 
174  // Calculate the extrapolated outlet patch flow rate
175  const scalar estimatedFlowRate = gSum(rhoOutlet*(patch().magSf()*nUp));
176 
177  if (estimatedFlowRate > 0.5*flowRate)
178  {
179  nUp *= (mag(flowRate)/mag(estimatedFlowRate));
180  }
181  else
182  {
183  nUp += ((flowRate - estimatedFlowRate)/gSum(rhoOutlet*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  operator==(Up);
191 }
192 
193 
195 {
196  if (updated())
197  {
198  return;
199  }
200 
201  // Find corresponding inlet patch
202  const label inletPatchID =
203  patch().patch().boundaryMesh().findPatchID(inletPatchName_);
204 
205  if (inletPatchID < 0)
206  {
208  << "Unable to find inlet patch " << inletPatchName_
209  << exit(FatalError);
210  }
211 
212  if (volumetric_)
213  {
214  updateValues(inletPatchID, one{}, one{});
215  }
216  else
217  {
218  // Mass flow-rate
219  if (db().foundObject<volScalarField>(rhoName_))
220  {
221  const volScalarField& rho = db().lookupObject<volScalarField>
222  (
223  rhoName_
224  );
225 
226  updateValues
227  (
228  inletPatchID,
229  rho.boundaryField()[patch().index()],
230  rho.boundaryField()[inletPatchID]
231  );
232  }
233  else
234  {
236  << "Cannot find density field " << rhoName_ << exit(FatalError);
237  }
238  }
239 
240  fixedValueFvPatchVectorField::updateCoeffs();
241 }
242 
243 
245 (
246  Ostream& os
247 ) const
248 {
250  os.writeEntry("inletPatch", inletPatchName_);
251  if (!volumetric_)
252  {
253  os.writeEntry("volumetric", volumetric_);
254  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
255  }
256  writeEntry("value", os);
257 }
258 
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 namespace Foam
263 {
265  (
268  );
269 }
270 
271 
272 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvPatch::Sf
const vectorField & Sf() const
Return face area vectors.
Definition: fvPatch.C:150
Foam::Ostream::writeEntryIfDifferent
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:248
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
matchedFlowRateOutletVelocityFvPatchVectorField.H
Foam::matchedFlowRateOutletVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: matchedFlowRateOutletVelocityFvPatchVectorField.C:245
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
rho
rho
Definition: readInitialConditions.H:88
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:80
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
one.H
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::matchedFlowRateOutletVelocityFvPatchVectorField::matchedFlowRateOutletVelocityFvPatchVectorField
matchedFlowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: matchedFlowRateOutletVelocityFvPatchVectorField.C:38
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::matchedFlowRateOutletVelocityFvPatchVectorField
Velocity outlet boundary condition which corrects the extrapolated velocity to match the flow rate of...
Definition: matchedFlowRateOutletVelocityFvPatchVectorField.H:95
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::matchedFlowRateOutletVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: matchedFlowRateOutletVelocityFvPatchVectorField.C:194
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::stringOps::evaluate
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.
Definition: stringOpsEvaluate.C:37