swirlFanVelocityFvPatchField.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) 2018-2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "unitConversion.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 void Foam::swirlFanVelocityFvPatchField::calcFanJump()
38 {
39  if (this->cyclicPatch().owner())
40  {
41  const scalar rpm = rpm_->value(this->db().time().timeOutputValue());
42 
43  const surfaceScalarField& phi =
44  db().lookupObject<surfaceScalarField>(phiName_);
45 
46  const fvPatchField<scalar>& pOwner =
47  patch().lookupPatchField<volScalarField, scalar>(pName_);
48 
49  const label nbrIndex = this->cyclicPatch().neighbPatchID();
50 
51  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrIndex];
52 
53  const fvPatchField<scalar>& pSlave =
54  nbrPatch.lookupPatchField<volScalarField, scalar>(pName_);
55 
56  scalarField deltaP(mag(pOwner - pSlave));
57 
58  if (phi.dimensions() == dimMass/dimTime)
59  {
60  deltaP /=
61  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
62  }
63 
64  const vector axisHat =
65  gSum(patch().nf()*patch().magSf())/gSum(patch().magSf());
66 
67  vectorField tanDir
68  (
69  axisHat ^ (patch().Cf() - origin_)
70  );
71 
72  tanDir /= (mag(tanDir) + SMALL);
73 
74  scalarField magTangU(patch().size(), Zero);
75 
76  if (useRealRadius_)
77  {
78  const vectorField& pCf = patch().Cf();
79 
80  forAll(pCf, i)
81  {
82  const scalar rMag = mag(pCf[i] - origin_);
83 
84  if (rMag > rInner_ && rMag < rOuter_)
85  {
86  magTangU[i] =
87  deltaP[i]/rMag/fanEff_/rpmToRads(rpm);
88  }
89  }
90  }
91  else
92  {
93  if (rEff_ <= 0)
94  {
96  << "Effective radius rEff should be specified in the "<< nl
97  << "dictionary." << nl
98  << exit(FatalError);
99  }
100  magTangU =
101  deltaP/rEff_/fanEff_/rpmToRads(rpm);
102  }
103 
104  // Calculate the tangential velocity
105  const vectorField tangentialVelocity(magTangU*tanDir);
106 
107  this->jump_ = tangentialVelocity;
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  const fvPatch& p,
118 )
119 :
121  phiName_("phi"),
122  pName_("p"),
123  rhoName_("rho"),
124  origin_(),
125  rpm_(),
126  rEff_(0.0),
127  fanEff_(1.0),
128  rInner_(0.0),
129  rOuter_(0.0),
130  useRealRadius_(false)
131 {}
132 
133 
135 (
136  const fvPatch& p,
138  const dictionary& dict
139 )
140 :
142  phiName_(dict.getOrDefault<word>("phi", "phi")),
143  pName_(dict.getOrDefault<word>("p", "p")),
144  rhoName_(dict.getOrDefault<word>("rho", "rho")),
145  origin_
146  (
148  (
149  "origin",
150  returnReduce(patch().size(), maxOp<label>())
151  ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
152  : Zero
153  )
154  ),
155  rpm_
156  (
157  this->cyclicPatch().owner()
158  ? Function1<scalar>::New("rpm", dict)
159  : nullptr
160  ),
161  rEff_(dict.getOrDefault<scalar>("rEff", 0)),
162  fanEff_(dict.getOrDefault<scalar>("fanEff", 1)),
163  rInner_(dict.getOrDefault<scalar>("rInner", 0)),
164  rOuter_(dict.getOrDefault<scalar>("rOuter", 0)),
165  useRealRadius_(dict.getOrDefault("useRealRadius", false))
166 {}
167 
168 
170 (
171  const swirlFanVelocityFvPatchField& ptf,
172  const fvPatch& p,
174  const fvPatchFieldMapper& mapper
175 )
176 :
177  fixedJumpFvPatchField<vector>(ptf, p, iF, mapper),
178  phiName_(ptf.phiName_),
179  pName_(ptf.pName_),
180  rhoName_(ptf.rhoName_),
181  origin_(ptf.origin_),
182  rpm_(ptf.rpm_.clone()),
183  rEff_(ptf.rEff_),
184  fanEff_(ptf.fanEff_),
185  rInner_(ptf.rInner_),
186  rOuter_(ptf.rOuter_),
187  useRealRadius_(ptf.useRealRadius_)
188 {}
189 
190 
192 (
194 )
195 :
197  phiName_(ptf.phiName_),
198  pName_(ptf.pName_),
199  rhoName_(ptf.rhoName_),
200  origin_(ptf.origin_),
201  rpm_(ptf.rpm_.clone()),
202  rEff_(ptf.rEff_),
203  rInner_(ptf.rInner_),
204  rOuter_(ptf.rOuter_),
205  useRealRadius_(ptf.useRealRadius_)
206 {}
207 
208 
210 (
211  const swirlFanVelocityFvPatchField& ptf,
213 )
214 :
216  phiName_(ptf.phiName_),
217  pName_(ptf.pName_),
218  rhoName_(ptf.rhoName_),
219  origin_(ptf.origin_),
220  rpm_(ptf.rpm_.clone()),
221  rEff_(ptf.rEff_),
222  rInner_(ptf.rInner_),
223  rOuter_(ptf.rOuter_),
224  useRealRadius_(ptf.useRealRadius_)
225 {}
226 
227 
228 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
229 
231 {
232  if (this->updated())
233  {
234  return;
235  }
236 
237  calcFanJump();
238 }
239 
240 
242 {
244 
245  if (this->cyclicPatch().owner())
246  {
247  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
248  os.writeEntryIfDifferent<word>("p", "p", pName_);
249  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
250  os.writeEntry("origin", origin_);
251 
252  if (rpm_)
253  {
254  rpm_->writeData(os);
255  }
256 
257  os.writeEntryIfDifferent<scalar>("rEff", 0.0, rEff_);
258  os.writeEntryIfDifferent<bool>("useRealRadius", false, useRealRadius_);
259  os.writeEntryIfDifferent<scalar>("rInner", 0.0, rInner_);
260  os.writeEntryIfDifferent<scalar>("rOuter", 0.0, rOuter_);
261  }
262 }
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 namespace Foam
268 {
270  (
273  );
274 }
275 
276 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
volFields.H
Foam::maxOp
Definition: ops.H:223
Foam::Ostream::writeEntryIfDifferent
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:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::swirlFanVelocityFvPatchField::swirlFanVelocityFvPatchField
swirlFanVelocityFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: swirlFanVelocityFvPatchField.C:115
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
unitConversion.H
Unit conversion functions.
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:86
swirlFanVelocityFvPatchField.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:203
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::fixedJumpFvPatchField< vector >
Foam::fixedJumpFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fixedJumpFvPatchField.C:157
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::cyclicFvPatchField< vector >::cyclicPatch
const cyclicFvPatch & cyclicPatch() const
Return local reference cast into the cyclic patch.
Definition: cyclicFvPatchField.H:165
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::swirlFanVelocityFvPatchField::write
virtual void write(Ostream &os) const
Write.
Definition: swirlFanVelocityFvPatchField.C:241
Foam::fvPatch::Cf
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:113
Foam::fvPatch::lookupPatchField
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Definition: fvPatchFvMeshTemplates.C:34
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::rpmToRads
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Definition: unitConversion.H:73
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::fvPatchField< vector >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::swirlFanVelocityFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: swirlFanVelocityFvPatchField.C:230
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::cyclicFvPatch::neighbPatchID
virtual label neighbPatchID() const
Return neighbour.
Definition: cyclicFvPatch.H:98
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField< vector >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::swirlFanVelocityFvPatchField
This boundary condition provides a jump condition for U across a cyclic pressure jump condition and a...
Definition: swirlFanVelocityFvPatchField.H:177
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fixedJumpFvPatchField< vector >::jump_
Field< vector > jump_
"jump" field
Definition: fixedJumpFvPatchField.H:108