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-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
37void 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 =
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
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 (
88 deltaP[i]
89 / stabilise
90 (
91 fanEff_ * rMag * rpmToRads(rpm),
92 VSMALL
93 )
94 );
95 }
96 }
97 }
98 else
99 {
100 if (rEff_ <= 0)
101 {
103 << "Effective radius 'rEff' was ill-specified in the "
104 << "dictionary." << nl
105 << exit(FatalError);
106 }
107 magTangU =
108 (
109 deltaP
110 / stabilise
111 (
112 fanEff_ * rEff_ * rpmToRads(rpm),
113 VSMALL
114 )
115 );
116 }
117
118 // Calculate the tangential velocity
119 const vectorField tangentialVelocity(magTangU*tanDir);
120
121 this->setJump(tangentialVelocity);
122 }
123}
124
125
126// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127
129(
130 const fvPatch& p,
132)
133:
135 phiName_("phi"),
136 pName_("p"),
137 rhoName_("rho"),
138 origin_(),
139 rpm_(nullptr),
140 fanEff_(1),
141 rEff_(0),
142 rInner_(0),
143 rOuter_(0),
144 useRealRadius_(false)
145{}
146
147
149(
150 const fvPatch& p,
152 const dictionary& dict
153)
154:
156 phiName_(dict.getOrDefault<word>("phi", "phi")),
157 pName_(dict.getOrDefault<word>("p", "p")),
158 rhoName_(dict.getOrDefault<word>("rho", "rho")),
159 origin_
160 (
161 dict.getOrDefault
162 (
163 "origin",
164 returnReduce(patch().size(), maxOp<label>())
165 ? gSum(patch().Cf()*patch().magSf())/gSum(patch().magSf())
166 : Zero
167 )
168 ),
169 rpm_
170 (
171 this->cyclicPatch().owner()
172 ? Function1<scalar>::New("rpm", dict, &db())
173 : nullptr
174 ),
175 fanEff_(dict.getOrDefault<scalar>("fanEff", 1)),
176 rEff_(dict.getOrDefault<scalar>("rEff", 0)),
177 rInner_(dict.getOrDefault<scalar>("rInner", 0)),
178 rOuter_(dict.getOrDefault<scalar>("rOuter", 0)),
179 useRealRadius_(dict.getOrDefault("useRealRadius", false))
180{}
181
182
184(
186 const fvPatch& p,
188 const fvPatchFieldMapper& mapper
189)
190:
191 fixedJumpFvPatchField<vector>(rhs, p, iF, mapper),
192 phiName_(rhs.phiName_),
193 pName_(rhs.pName_),
194 rhoName_(rhs.rhoName_),
195 origin_(rhs.origin_),
196 rpm_(rhs.rpm_.clone()),
197 fanEff_(rhs.fanEff_),
198 rEff_(rhs.rEff_),
199 rInner_(rhs.rInner_),
200 rOuter_(rhs.rOuter_),
201 useRealRadius_(rhs.useRealRadius_)
202{}
203
204
206(
208)
209:
211 phiName_(rhs.phiName_),
212 pName_(rhs.pName_),
213 rhoName_(rhs.rhoName_),
214 origin_(rhs.origin_),
215 rpm_(rhs.rpm_.clone()),
216 fanEff_(rhs.fanEff_),
217 rEff_(rhs.rEff_),
218 rInner_(rhs.rInner_),
219 rOuter_(rhs.rOuter_),
220 useRealRadius_(rhs.useRealRadius_)
221{}
222
223
225(
228)
229:
231 phiName_(rhs.phiName_),
232 pName_(rhs.pName_),
233 rhoName_(rhs.rhoName_),
234 origin_(rhs.origin_),
235 rpm_(rhs.rpm_.clone()),
236 fanEff_(rhs.fanEff_),
237 rEff_(rhs.rEff_),
238 rInner_(rhs.rInner_),
239 rOuter_(rhs.rOuter_),
240 useRealRadius_(rhs.useRealRadius_)
241{}
242
243
244// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
245
247{
248 if (this->updated())
249 {
250 return;
251 }
252
253 calcFanJump();
254}
255
256
258{
260
261 if (this->cyclicPatch().owner())
262 {
263 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
264 os.writeEntryIfDifferent<word>("p", "p", pName_);
265 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
266 os.writeEntry("origin", origin_);
267
268 if (rpm_)
269 {
270 rpm_->writeData(os);
271 }
272
273 os.writeEntryIfDifferent<scalar>("fanEff", 1, fanEff_);
274
275 if (useRealRadius_)
276 {
277 os.writeEntry("useRealRadius", "true");
278 os.writeEntryIfDifferent<scalar>("rInner", 0, rInner_);
279 os.writeEntryIfDifferent<scalar>("rOuter", 0, rOuter_);
280 }
281 else
282 {
283 os.writeEntryIfDifferent<scalar>("rEff", 0, rEff_);
284 }
285 }
286}
287
288
289// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290
291namespace Foam
292{
294 (
297 );
298}
299
300// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
label size() const noexcept
The number of elements in the UList.
Definition: UListI.H:420
const cyclicFvPatch & cyclicPatch() const
Return local reference cast into the cyclic patch.
virtual label neighbPatchID() const
Return neighbour.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void setJump(const Field< vector > &jump)
Set the jump field.
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:209
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:119
const Type & lookupObject(const word &name, const bool recursive=false) const
This boundary condition provides a jump condition for U across a cyclic pressure jump condition and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar rpmToRads() noexcept
Multiplication factor for revolutions/minute to radians/sec.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.
Unit conversion functions.