swirlFlowRateInletVelocityFvPatchVectorField.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) 2011-2021 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"
29 #include "fvPatchFieldMapper.H"
30 #include "surfaceFields.H"
31 #include "mathematicalConstants.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
43  phiName_("phi"),
44  rhoName_("rho"),
45  origin_(),
46  axis_(Zero),
47  flowRate_(),
48  rpm_()
49 {}
50 
51 
54 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
61  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
62  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
63  origin_
64  (
65  dict.lookupOrDefault
66  (
67  "origin",
68  returnReduce(patch().size(), sumOp<label>())
69  ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
70  : Zero
71  )
72  ),
73  axis_
74  (
75  dict.lookupOrDefault
76  (
77  "axis",
78  returnReduce(patch().size(), sumOp<label>())
79  ? -gSum(patch().Sf())/gSum(patch().magSf())
80  : Zero
81  )
82  ),
83  flowRate_(Function1<scalar>::New("flowRate", dict)),
84  rpm_(Function1<scalar>::New("rpm", dict))
85 {}
86 
87 
90 (
92  const fvPatch& p,
94  const fvPatchFieldMapper& mapper
95 )
96 :
97  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
98  phiName_(ptf.phiName_),
99  rhoName_(ptf.rhoName_),
100  origin_(ptf.origin_),
101  axis_(ptf.axis_),
102  flowRate_(ptf.flowRate_, false),
103  rpm_(ptf.rpm_, false)
104 {}
105 
106 
109 (
112 )
113 :
115  phiName_(ptf.phiName_),
116  rhoName_(ptf.rhoName_),
117  origin_(ptf.origin_),
118  axis_(ptf.axis_),
119  flowRate_(ptf.flowRate_, false),
120  rpm_(ptf.rpm_, false)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 {
128  if (updated())
129  {
130  return;
131  }
132 
133  const scalar t = this->db().time().userTimeValue();
134  const scalar flowRate = flowRate_->value(t);
135  const scalar rpm = rpm_->value(t);
136 
137  const scalar totArea = gSum(patch().magSf());
138  const scalar avgU = -flowRate/totArea;
139 
140  const vector axisHat = axis_/mag(axis_);
141 
142  // Update angular velocity - convert [rpm] to [rad/s]
143  tmp<vectorField> tangentialVelocity
144  (
145  axisHat ^ (rpm*constant::mathematical::pi/30.0)*(patch().Cf() - origin_)
146  );
147 
148  tmp<vectorField> n = patch().nf();
149 
150  const surfaceScalarField& phi =
151  db().lookupObject<surfaceScalarField>(phiName_);
152 
153  if (phi.dimensions() == dimFlux)
154  {
155  // volumetric flow-rate
156  operator==(tangentialVelocity + n*avgU);
157  }
158  else if (phi.dimensions() == dimMassFlux)
159  {
160  const fvPatchField<scalar>& rhop =
161  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
162 
163  // mass flow-rate
164  operator==(tangentialVelocity + n*avgU/rhop);
165  }
166  else
167  {
169  << "dimensions of " << phiName_ << " are incorrect" << nl
170  << " on patch " << this->patch().name()
171  << " of field " << this->internalField().name()
172  << " in file " << this->internalField().objectPath()
173  << nl << exit(FatalError);
174  }
175 
177 }
178 
179 
181 (
182  Ostream& os
183 ) const
184 {
186  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
187  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
188  writeEntry(os, "origin", origin_);
189  writeEntry(os, "axis", axis_);
190  writeEntry(os, flowRate_());
191  writeEntry(os, rpm_());
192  writeEntry(os, "value", *this);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 namespace Foam
199 {
201  (
204  );
205 }
206 
207 
208 // ************************************************************************* //
Foam::surfaceFields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Run-time selectable general function of one variable.
Definition: Function1.H:52
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet & dimensions() const
Return dimensions.
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
swirlFlowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
phi
Definition: correctPhi.H:3
static const zero Zero
Definition: zero.H:97
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition provides a volumetric- OR mass-flow normal vector boundary condition by its m...
Namespace for OpenFOAM.