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-2018 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  patch().size()
69  ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
70  : Zero
71  )
72  ),
73  axis_
74  (
75  dict.lookupOrDefault
76  (
77  "axis",
78  patch().size()
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 (
111 )
112 :
114  phiName_(ptf.phiName_),
115  rhoName_(ptf.rhoName_),
116  origin_(ptf.origin_),
117  axis_(ptf.axis_),
118  flowRate_(ptf.flowRate_, false),
119  rpm_(ptf.rpm_, false)
120 {}
121 
122 
125 (
128 )
129 :
131  phiName_(ptf.phiName_),
132  rhoName_(ptf.rhoName_),
133  origin_(ptf.origin_),
134  axis_(ptf.axis_),
135  flowRate_(ptf.flowRate_, false),
136  rpm_(ptf.rpm_, false)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  const scalar t = this->db().time().timeOutputValue();
150  const scalar flowRate = flowRate_->value(t);
151  const scalar rpm = rpm_->value(t);
152 
153  const scalar totArea = gSum(patch().magSf());
154  const scalar avgU = -flowRate/totArea;
155 
156  const vector axisHat = axis_/mag(axis_);
157 
158  // Update angular velocity - convert [rpm] to [rad/s]
159  tmp<vectorField> tangentialVelocity
160  (
161  axisHat ^ (rpm*constant::mathematical::pi/30.0)*(patch().Cf() - origin_)
162  );
163 
164  tmp<vectorField> n = patch().nf();
165 
166  const surfaceScalarField& phi =
167  db().lookupObject<surfaceScalarField>(phiName_);
168 
169  if (phi.dimensions() == dimVelocity*dimArea)
170  {
171  // volumetric flow-rate
172  operator==(tangentialVelocity + n*avgU);
173  }
174  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
175  {
176  const fvPatchField<scalar>& rhop =
177  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
178 
179  // mass flow-rate
180  operator==(tangentialVelocity + n*avgU/rhop);
181  }
182  else
183  {
185  << "dimensions of " << phiName_ << " are incorrect" << nl
186  << " on patch " << this->patch().name()
187  << " of field " << this->internalField().name()
188  << " in file " << this->internalField().objectPath()
189  << nl << exit(FatalError);
190  }
191 
193 }
194 
195 
197 (
198  Ostream& os
199 ) const
200 {
202  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
203  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
204  os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
205  os.writeKeyword("axis") << axis_ << token::END_STATEMENT << nl;
206  flowRate_->writeData(os);
207  rpm_->writeData(os);
208  writeEntry("value", os);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 namespace Foam
215 {
217  (
220  );
221 }
222 
223 
224 // ************************************************************************* //
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
surfaceScalarField & phi
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:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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:362
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
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.
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:53
static const char nl
Definition: Ostream.H:265
const dimensionSet dimDensity
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:311
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
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...
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet dimVelocity