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-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"
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  returnReduce(patch().size(), maxOp<label>())
72  ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
73  : Zero
74  )
75  ),
76  axis_
77  (
78  dict.getOrDefault
79  (
80  "axis",
81  returnReduce(patch().size(), maxOp<label>())
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 :
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 surfaceScalarField& phi =
172  db().lookupObject<surfaceScalarField>(phiName_);
173 
174  if (phi.dimensions() == dimVelocity*dimArea)
175  {
176  // volumetric flow-rate
177  operator==(tangentialVelocity + n*avgU);
178  }
179  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
180  {
181  const fvPatchField<scalar>& rhop =
182  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
183 
184  // mass flow-rate
185  operator==(tangentialVelocity + n*avgU/rhop);
186  }
187  else
188  {
190  << "dimensions of " << phiName_ << " are incorrect" << nl
191  << " on patch " << this->patch().name()
192  << " of field " << this->internalField().name()
193  << " in file " << this->internalField().objectPath()
194  << nl << exit(FatalError);
195  }
196  }
197 
199 }
200 
201 
203 (
204  Ostream& os
205 ) const
206 {
208  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
209  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
210  os.writeEntry("origin", origin_);
211  os.writeEntry("axis", axis_);
212  flowRate_->writeData(os);
213  rpm_->writeData(os);
214  writeEntry("value", os);
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 namespace Foam
221 {
223  (
226  );
227 }
228 
229 
230 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::maxOp
Definition: ops.H:223
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
Foam::swirlFlowRateInletVelocityFvPatchVectorField::swirlFlowRateInletVelocityFvPatchVectorField
swirlFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: swirlFlowRateInletVelocityFvPatchVectorField.C:40
Foam::swirlFlowRateInletVelocityFvPatchVectorField
This boundary condition provides a volumetric- OR mass-flow normal vector boundary condition by its m...
Definition: swirlFlowRateInletVelocityFvPatchVectorField.H:126
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::swirlFlowRateInletVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: swirlFlowRateInletVelocityFvPatchVectorField.C:203
unitConversion.H
Unit conversion functions.
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
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::swirlFlowRateInletVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: swirlFlowRateInletVelocityFvPatchVectorField.C:145
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
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
swirlFlowRateInletVelocityFvPatchVectorField.H
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::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
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
Foam::rpmToRads
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Definition: unitConversion.H:73
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
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::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54