swirlFlowRateInletVelocityFvPatchVectorField.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2022 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"
32 #include "fvPatchFieldMapper.H"
33 #include "surfaceFields.H"
34 #include "unitConversion.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
46  phiName_("phi"),
47  rhoName_("rho"),
48  origin_(),
49  axis_(Zero),
50  flowRate_(),
51  rpm_()
52 {}
53 
54 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
64  phiName_(dict.getOrDefault<word>("phi", "phi")),
65  rhoName_(dict.getOrDefault<word>("rho", "rho")),
66  origin_
67  (
68  dict.getOrDefault
69  (
70  "origin",
71  returnReduceOr(patch().size())
72  ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
73  : Zero
74  )
75  ),
76  axis_
77  (
78  dict.getOrDefault
79  (
80  "axis",
81  returnReduceOr(patch().size())
82  ? -gSum(patch().Sf())/gSum(patch().magSf())
83  : Zero
84  )
85  ),
86  flowRate_(Function1<scalar>::New("flowRate", dict, &db())),
87  rpm_(Function1<scalar>::New("rpm", dict, &db()))
88 {}
89 
90 
93 (
95  const fvPatch& p,
97  const fvPatchFieldMapper& mapper
98 )
99 :
100  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
101  phiName_(ptf.phiName_),
102  rhoName_(ptf.rhoName_),
103  origin_(ptf.origin_),
104  axis_(ptf.axis_),
105  flowRate_(ptf.flowRate_.clone()),
106  rpm_(ptf.rpm_.clone())
107 {}
108 
109 
112 (
114 )
115 :
117  phiName_(ptf.phiName_),
118  rhoName_(ptf.rhoName_),
119  origin_(ptf.origin_),
120  axis_(ptf.axis_),
121  flowRate_(ptf.flowRate_.clone()),
122  rpm_(ptf.rpm_.clone())
123 {}
124 
125 
128 (
131 )
132 :
133  fixedValueFvPatchField<vector>(ptf, iF),
134  phiName_(ptf.phiName_),
135  rhoName_(ptf.rhoName_),
136  origin_(ptf.origin_),
137  axis_(ptf.axis_),
138  flowRate_(ptf.flowRate_.clone()),
139  rpm_(ptf.rpm_.clone())
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
146 {
147  if (updated())
148  {
149  return;
150  }
151  const scalar totArea = gSum(patch().magSf());
152 
153  if (totArea > ROOTVSMALL && axis_ != vector(Zero))
154  {
155  const scalar t = this->db().time().timeOutputValue();
156  const scalar flowRate = flowRate_->value(t);
157  const scalar omega = rpmToRads(rpm_->value(t));
158 
159  const scalar avgU = -flowRate/totArea;
160 
161  const vector axisHat = axis_/mag(axis_);
162 
163  // Update angular velocity
164  tmp<vectorField> tangentialVelocity
165  (
166  axisHat ^ omega*(patch().Cf() - origin_)
167  );
168 
169  tmp<vectorField> n = patch().nf();
170 
171  const auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
172 
173  if (phi.dimensions() == dimVolume/dimTime)
174  {
175  // volumetric flow-rate
176  operator==(tangentialVelocity + n*avgU);
177  }
178  else if (phi.dimensions() == dimMass/dimTime)
179  {
180  const auto& rhop =
181  patch().lookupPatchField<volScalarField>(rhoName_);
182 
183  // mass flow-rate
184  operator==(tangentialVelocity + n*avgU/rhop);
185  }
186  else
187  {
189  << "dimensions of " << phiName_ << " are incorrect" << nl
190  << " on patch " << this->patch().name()
191  << " of field " << this->internalField().name()
192  << " in file " << this->internalField().objectPath()
193  << nl << exit(FatalError);
194  }
195  }
196 
198 }
199 
200 
202 (
203  Ostream& os
204 ) const
205 {
207  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
208  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
209  os.writeEntry("origin", origin_);
210  os.writeEntry("axis", axis_);
211  flowRate_->writeData(os);
212  rpm_->writeData(os);
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 namespace Foam
220 {
222  (
224  swirlFlowRateInletVelocityFvPatchVectorField
225  );
226 }
227 
228 
229 // ************************************************************************* //
Foam::surfaceFields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
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
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
swirlFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A FieldMapper for finite-volume patch fields.
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)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
label n
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition provides a volumetric- OR mass-flow normal vector boundary condition by its m...
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127