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-------------------------------------------------------------------------------
11License
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
38Foam::swirlFlowRateInletVelocityFvPatchVectorField::
39swirlFlowRateInletVelocityFvPatchVectorField
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
55Foam::swirlFlowRateInletVelocityFvPatchVectorField::
56swirlFlowRateInletVelocityFvPatchVectorField
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
91Foam::swirlFlowRateInletVelocityFvPatchVectorField::
92swirlFlowRateInletVelocityFvPatchVectorField
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
110Foam::swirlFlowRateInletVelocityFvPatchVectorField::
111swirlFlowRateInletVelocityFvPatchVectorField
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
126Foam::swirlFlowRateInletVelocityFvPatchVectorField::
127swirlFlowRateInletVelocityFvPatchVectorField
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
175 {
176 // volumetric flow-rate
177 operator==(tangentialVelocity + n*avgU);
178 }
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
220namespace Foam
221{
223 (
226 );
227}
228
229
230// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a volumetric- OR mass-flow normal vector boundary condition by its m...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
constexpr scalar rpmToRads() noexcept
Multiplication factor for revolutions/minute to radians/sec.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
Vector< scalar > vector
Definition: vector.H:61
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
Foam::surfaceFields.
Unit conversion functions.