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 | Copyright (C) 2011-2017 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  gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
69  )
70  ),
71  axis_
72  (
73  dict.lookupOrDefault
74  (
75  "axis",
76  -gSum(patch().Sf())/gSum(patch().magSf())
77  )
78  ),
79  flowRate_(Function1<scalar>::New("flowRate", dict)),
80  rpm_(Function1<scalar>::New("rpm", dict))
81 {}
82 
83 
86 (
88  const fvPatch& p,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
94  phiName_(ptf.phiName_),
95  rhoName_(ptf.rhoName_),
96  origin_(ptf.origin_),
97  axis_(ptf.axis_),
98  flowRate_(ptf.flowRate_, false),
99  rpm_(ptf.rpm_, false)
100 {}
101 
102 
105 (
107 )
108 :
110  phiName_(ptf.phiName_),
111  rhoName_(ptf.rhoName_),
112  origin_(ptf.origin_),
113  axis_(ptf.axis_),
114  flowRate_(ptf.flowRate_, false),
115  rpm_(ptf.rpm_, false)
116 {}
117 
118 
121 (
124 )
125 :
127  phiName_(ptf.phiName_),
128  rhoName_(ptf.rhoName_),
129  origin_(ptf.origin_),
130  axis_(ptf.axis_),
131  flowRate_(ptf.flowRate_, false),
132  rpm_(ptf.rpm_, false)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  if (updated())
141  {
142  return;
143  }
144 
145  const scalar t = this->db().time().timeOutputValue();
146  const scalar flowRate = flowRate_->value(t);
147  const scalar rpm = rpm_->value(t);
148 
149  const scalar totArea = gSum(patch().magSf());
150  const scalar avgU = -flowRate/totArea;
151 
152  const vector axisHat = axis_/mag(axis_);
153 
154  // Update angular velocity - convert [rpm] to [rad/s]
155  tmp<vectorField> tangentialVelocity
156  (
157  axisHat ^ (rpm*constant::mathematical::pi/30.0)*(patch().Cf() - origin_)
158  );
159 
160  tmp<vectorField> n = patch().nf();
161 
162  const surfaceScalarField& phi =
163  db().lookupObject<surfaceScalarField>(phiName_);
164 
165  if (phi.dimensions() == dimVelocity*dimArea)
166  {
167  // volumetric flow-rate
168  operator==(tangentialVelocity + n*avgU);
169  }
170  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
171  {
172  const fvPatchField<scalar>& rhop =
173  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
174 
175  // mass flow-rate
176  operator==(tangentialVelocity + n*avgU/rhop);
177  }
178  else
179  {
181  << "dimensions of " << phiName_ << " are incorrect" << nl
182  << " on patch " << this->patch().name()
183  << " of field " << this->internalField().name()
184  << " in file " << this->internalField().objectPath()
185  << nl << exit(FatalError);
186  }
187 
189 }
190 
191 
193 (
194  Ostream& os
195 ) const
196 {
198  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
199  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
200  os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
201  os.writeKeyword("axis") << axis_ << token::END_STATEMENT << nl;
202  flowRate_->writeData(os);
203  rpm_->writeData(os);
204  writeEntry("value", os);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 namespace Foam
211 {
213  (
216  );
217 }
218 
219 
220 // ************************************************************************* //
Foam::surfaceFields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
surfaceScalarField & phi
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:53
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: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
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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:91
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:262
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
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