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-2022 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 auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
44 
45  const auto& pOwner = patch().lookupPatchField<volScalarField>(pName_);
46 
47  const label nbrIndex = this->cyclicPatch().neighbPatchID();
48 
49  const fvPatch& nbrPatch = patch().boundaryMesh()[nbrIndex];
50 
51  const auto& pSlave = nbrPatch.lookupPatchField<volScalarField>(pName_);
52 
53  scalarField deltaP(mag(pOwner - pSlave));
54 
55  if (phi.dimensions() == dimMass/dimTime)
56  {
57  deltaP /= patch().lookupPatchField<volScalarField>(rhoName_);
58  }
59 
60  const vector axisHat =
61  gSum(patch().nf()*patch().magSf())/gSum(patch().magSf());
62 
63  vectorField tanDir
64  (
65  axisHat ^ (patch().Cf() - origin_)
66  );
67 
68  tanDir /= (mag(tanDir) + SMALL);
69 
70  scalarField magTangU(patch().size(), Zero);
71 
72  if (useRealRadius_)
73  {
74  const vectorField& pCf = patch().Cf();
75 
76  forAll(pCf, i)
77  {
78  const scalar rMag = mag(pCf[i] - origin_);
79 
80  if (rMag > rInner_ && rMag < rOuter_)
81  {
82  magTangU[i] =
83  (
84  deltaP[i]
85  / stabilise
86  (
87  fanEff_ * rMag * rpmToRads(rpm),
88  VSMALL
89  )
90  );
91  }
92  }
93  }
94  else
95  {
96  if (rEff_ <= 0)
97  {
99  << "Effective radius 'rEff' was ill-specified in the "
100  << "dictionary." << nl
101  << exit(FatalError);
102  }
103  magTangU =
104  (
105  deltaP
106  / stabilise
107  (
108  fanEff_ * rEff_ * rpmToRads(rpm),
109  VSMALL
110  )
111  );
112  }
113 
114  // Calculate the tangential velocity
115  const vectorField tangentialVelocity(magTangU*tanDir);
116 
117  this->setJump(tangentialVelocity);
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
125 (
126  const fvPatch& p,
127  const DimensionedField<vector, volMesh>& iF
128 )
129 :
130  fixedJumpFvPatchField<vector>(p, iF),
131  phiName_("phi"),
132  pName_("p"),
133  rhoName_("rho"),
134  origin_(),
135  rpm_(nullptr),
136  fanEff_(1),
137  rEff_(0),
138  rInner_(0),
139  rOuter_(0),
140  useRealRadius_(false)
141 {}
142 
143 
145 (
146  const fvPatch& p,
148  const dictionary& dict
149 )
150 :
152  phiName_(dict.getOrDefault<word>("phi", "phi")),
153  pName_(dict.getOrDefault<word>("p", "p")),
154  rhoName_(dict.getOrDefault<word>("rho", "rho")),
155  origin_
156  (
157  dict.getOrDefault
158  (
159  "origin",
160  returnReduceOr(patch().size())
161  ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
162  : Zero
163  )
164  ),
165  rpm_
166  (
167  this->cyclicPatch().owner()
168  ? Function1<scalar>::New("rpm", dict, &db())
169  : nullptr
170  ),
171  fanEff_(dict.getOrDefault<scalar>("fanEff", 1)),
172  rEff_(dict.getOrDefault<scalar>("rEff", 0)),
173  rInner_(dict.getOrDefault<scalar>("rInner", 0)),
174  rOuter_(dict.getOrDefault<scalar>("rOuter", 0)),
175  useRealRadius_(dict.getOrDefault("useRealRadius", false))
176 {}
177 
178 
180 (
181  const swirlFanVelocityFvPatchField& rhs,
182  const fvPatch& p,
184  const fvPatchFieldMapper& mapper
185 )
186 :
187  fixedJumpFvPatchField<vector>(rhs, p, iF, mapper),
188  phiName_(rhs.phiName_),
189  pName_(rhs.pName_),
190  rhoName_(rhs.rhoName_),
191  origin_(rhs.origin_),
192  rpm_(rhs.rpm_.clone()),
193  fanEff_(rhs.fanEff_),
194  rEff_(rhs.rEff_),
195  rInner_(rhs.rInner_),
196  rOuter_(rhs.rOuter_),
197  useRealRadius_(rhs.useRealRadius_)
198 {}
199 
200 
202 (
204 )
205 :
207  phiName_(rhs.phiName_),
208  pName_(rhs.pName_),
209  rhoName_(rhs.rhoName_),
210  origin_(rhs.origin_),
211  rpm_(rhs.rpm_.clone()),
212  fanEff_(rhs.fanEff_),
213  rEff_(rhs.rEff_),
214  rInner_(rhs.rInner_),
215  rOuter_(rhs.rOuter_),
216  useRealRadius_(rhs.useRealRadius_)
217 {}
218 
219 
221 (
222  const swirlFanVelocityFvPatchField& rhs,
224 )
225 :
226  fixedJumpFvPatchField<vector>(rhs, iF),
227  phiName_(rhs.phiName_),
228  pName_(rhs.pName_),
229  rhoName_(rhs.rhoName_),
230  origin_(rhs.origin_),
231  rpm_(rhs.rpm_.clone()),
232  fanEff_(rhs.fanEff_),
233  rEff_(rhs.rEff_),
234  rInner_(rhs.rInner_),
235  rOuter_(rhs.rOuter_),
236  useRealRadius_(rhs.useRealRadius_)
237 {}
238 
239 
240 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241 
243 {
244  if (this->updated())
245  {
246  return;
247  }
248 
249  calcFanJump();
250 }
251 
252 
254 {
256 
257  if (this->cyclicPatch().owner())
258  {
259  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
260  os.writeEntryIfDifferent<word>("p", "p", pName_);
261  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
262  os.writeEntry("origin", origin_);
263 
264  if (rpm_)
265  {
266  rpm_->writeData(os);
267  }
268 
269  os.writeEntryIfDifferent<scalar>("fanEff", 1, fanEff_);
270 
271  if (useRealRadius_)
272  {
273  os.writeEntry("useRealRadius", "true");
274  os.writeEntryIfDifferent<scalar>("rInner", 0, rInner_);
275  os.writeEntryIfDifferent<scalar>("rOuter", 0, rOuter_);
276  }
277  else
278  {
279  os.writeEntryIfDifferent<scalar>("rEff", 0, rEff_);
280  }
281  }
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 namespace Foam
288 {
290  (
292  swirlFanVelocityFvPatchField
293  );
294 }
295 
296 // ************************************************************************* //
Foam::surfaceFields.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
This boundary condition provides a jump condition for U across a cyclic pressure jump condition and a...
virtual label neighbPatchID() const
Return neighbour.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const objectRegistry & db() const
The associated objectRegistry.
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr scalar rpmToRads(const scalar rpm) noexcept
Conversion from revolutions/minute to radians/sec.
Unit conversion functions.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const AnyType *=nullptr) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Type gSum(const FieldField< Field, Type > &f)
const cyclicFvPatch & cyclicPatch() const
Return local reference cast into the cyclic patch.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:113
Vector< scalar > vector
Definition: vector.H:57
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:336
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition: fvPatch.H:260
swirlFanVelocityFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void write(Ostream &) const
Write.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
OBJstream os(runTime.globalPath()/outputName)
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
This boundary condition provides a jump condition, using the cyclic condition as a base...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual void setJump(const Field< vector > &jump)
Set the jump field.
Field< vector > vectorField
Specialisation of Field<T> for vector.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127