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-2016 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  flowRate_(),
46  rpm_()
47 {}
48 
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
60  phiName_(ptf.phiName_),
61  rhoName_(ptf.rhoName_),
62  flowRate_(ptf.flowRate_, false),
63  rpm_(ptf.rpm_, false)
64 {}
65 
66 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
76  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
77  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
78  flowRate_(Function1<scalar>::New("flowRate", dict)),
79  rpm_(Function1<scalar>::New("rpm", dict))
80 {}
81 
82 
85 (
87 )
88 :
90  phiName_(ptf.phiName_),
91  rhoName_(ptf.rhoName_),
92  flowRate_(ptf.flowRate_, false),
93  rpm_(ptf.rpm_, false)
94 {}
95 
96 
99 (
102 )
103 :
105  phiName_(ptf.phiName_),
106  rhoName_(ptf.rhoName_),
107  flowRate_(ptf.flowRate_, false),
108  rpm_(ptf.rpm_, false)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  if (updated())
117  {
118  return;
119  }
120 
121  const scalar t = this->db().time().timeOutputValue();
122  const scalar flowRate = flowRate_->value(t);
123  const scalar rpm = rpm_->value(t);
124 
125  const scalar totArea = gSum(patch().magSf());
126  const scalar avgU = -flowRate/totArea;
127 
128  const vector avgCenter = gSum(patch().Cf()*patch().magSf())/totArea;
129  const vector avgNormal = gSum(patch().Sf())/totArea;
130 
131  // Update angular velocity - convert [rpm] to [rad/s]
132  tmp<vectorField> tangentialVelocity
133  (
134  (rpm*constant::mathematical::pi/30.0)
135  * (patch().Cf() - avgCenter) ^ avgNormal
136  );
137 
138  tmp<vectorField> n = patch().nf();
139 
140  const surfaceScalarField& phi =
141  db().lookupObject<surfaceScalarField>(phiName_);
142 
143  if (phi.dimensions() == dimVelocity*dimArea)
144  {
145  // volumetric flow-rate
146  operator==(tangentialVelocity + n*avgU);
147  }
148  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
149  {
150  const fvPatchField<scalar>& rhop =
151  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
152 
153  // mass flow-rate
154  operator==(tangentialVelocity + n*avgU/rhop);
155  }
156  else
157  {
159  << "dimensions of " << phiName_ << " are incorrect" << nl
160  << " on patch " << this->patch().name()
161  << " of field " << this->internalField().name()
162  << " in file " << this->internalField().objectPath()
163  << nl << exit(FatalError);
164  }
165 
167 }
168 
169 
171 (
172  Ostream& os
173 ) const
174 {
176  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
177  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
178  flowRate_->writeData(os);
179  rpm_->writeData(os);
180  writeEntry("value", os);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 namespace Foam
187 {
189  (
192  );
193 }
194 
195 
196 // ************************************************************************* //
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:65
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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 & dimensions() const
Return dimensions.
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:312
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
label n
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:54
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
const dimensionSet dimVelocity