Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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 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 = gSum(patch().Sf())/gSum(patch().magSf());
61
62 vectorField tanDir
63 (
64 axisHat ^ (patch().Cf() - origin_)
65 );
66
67 tanDir /= (mag(tanDir) + SMALL);
68
69 scalarField magTangU(patch().size(), Zero);
70
71 if (useRealRadius_)
72 {
73 const vectorField& pCf = patch().Cf();
74
75 forAll(pCf, i)
76 {
77 const scalar rMag = mag(pCf[i] - origin_);
78
79 if (rMag > rInner_ && rMag < rOuter_)
80 {
81 magTangU[i] =
82 (
83 deltaP[i]
84 / stabilise
85 (
86 fanEff_ * rMag * rpmToRads(rpm),
87 VSMALL
88 )
89 );
90 }
91 }
92 }
93 else
94 {
95 if (rEff_ <= 0)
96 {
98 << "Effective radius 'rEff' was ill-specified in the "
99 << "dictionary." << nl
100 << exit(FatalError);
101 }
102 magTangU =
103 (
104 deltaP
105 / stabilise
106 (
107 fanEff_ * rEff_ * rpmToRads(rpm),
108 VSMALL
109 )
110 );
111 }
112
113 // Calculate the tangential velocity
114 const vectorField tangentialVelocity(magTangU*tanDir);
115
116 this->setJump(tangentialVelocity);
117 }
118}
119
120
121// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122
124(
125 const fvPatch& p,
126 const DimensionedField<vector, volMesh>& iF
127)
128:
129 fixedJumpFvPatchField<vector>(p, iF),
130 phiName_("phi"),
131 pName_("p"),
132 rhoName_("rho"),
133 origin_(Zero),
134 rpm_(nullptr),
135 fanEff_(1),
136 rEff_(0),
137 rInner_(0),
138 rOuter_(0),
139 useRealRadius_(false)
140{}
141
142
144(
145 const fvPatch& p,
147 const dictionary& dict
148)
149:
151 phiName_(dict.getOrDefault<word>("phi", "phi")),
152 pName_(dict.getOrDefault<word>("p", "p")),
153 rhoName_(dict.getOrDefault<word>("rho", "rho")),
154 origin_
155 (
156 dict.getOrDefault
157 (
158 "origin",
159 returnReduceOr(patch().size())
160 ? gWeightedAverage(patch().magSf(), patch().Cf())
161 : Zero
162 )
163 ),
164 rpm_
165 (
166 this->cyclicPatch().owner()
167 ? Function1<scalar>::New("rpm", dict, &db())
168 : nullptr
169 ),
170 fanEff_(dict.getOrDefault<scalar>("fanEff", 1)),
171 rEff_(dict.getOrDefault<scalar>("rEff", 0)),
172 rInner_(dict.getOrDefault<scalar>("rInner", 0)),
173 rOuter_(dict.getOrDefault<scalar>("rOuter", 0)),
174 useRealRadius_(dict.getOrDefault("useRealRadius", false))
175{}
176
177
179(
181 const fvPatch& p,
183 const fvPatchFieldMapper& mapper
184)
185:
186 fixedJumpFvPatchField<vector>(rhs, p, iF, mapper),
187 phiName_(rhs.phiName_),
188 pName_(rhs.pName_),
189 rhoName_(rhs.rhoName_),
190 origin_(rhs.origin_),
191 rpm_(rhs.rpm_.clone()),
192 fanEff_(rhs.fanEff_),
193 rEff_(rhs.rEff_),
194 rInner_(rhs.rInner_),
195 rOuter_(rhs.rOuter_),
196 useRealRadius_(rhs.useRealRadius_)
197{}
198
199
201(
203)
204:
206 phiName_(rhs.phiName_),
207 pName_(rhs.pName_),
208 rhoName_(rhs.rhoName_),
209 origin_(rhs.origin_),
210 rpm_(rhs.rpm_.clone()),
211 fanEff_(rhs.fanEff_),
212 rEff_(rhs.rEff_),
213 rInner_(rhs.rInner_),
214 rOuter_(rhs.rOuter_),
215 useRealRadius_(rhs.useRealRadius_)
216{}
217
218
220(
223)
224:
226 phiName_(rhs.phiName_),
227 pName_(rhs.pName_),
228 rhoName_(rhs.rhoName_),
229 origin_(rhs.origin_),
230 rpm_(rhs.rpm_.clone()),
231 fanEff_(rhs.fanEff_),
232 rEff_(rhs.rEff_),
233 rInner_(rhs.rInner_),
234 rOuter_(rhs.rOuter_),
235 useRealRadius_(rhs.useRealRadius_)
236{}
237
238
239// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
240
242{
243 if (this->updated())
244 {
245 return;
246 }
247
248 calcFanJump();
249}
250
251
253{
255
256 if (this->cyclicPatch().owner())
257 {
258 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
259 os.writeEntryIfDifferent<word>("p", "p", pName_);
260 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
261 os.writeEntry("origin", origin_);
262
263 if (rpm_)
264 {
265 rpm_->writeData(os);
266 }
267
268 os.writeEntryIfDifferent<scalar>("fanEff", 1, fanEff_);
269
270 if (useRealRadius_)
271 {
272 os.writeEntry("useRealRadius", "true");
273 os.writeEntryIfDifferent<scalar>("rInner", 0, rInner_);
274 os.writeEntryIfDifferent<scalar>("rOuter", 0, rOuter_);
275 }
276 else
277 {
278 os.writeEntryIfDifferent<scalar>("rEff", 0, rEff_);
279 }
280 }
282
283
284// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285
286namespace Foam
287{
289 (
291 swirlFanVelocityFvPatchField
292 );
293}
294
295// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition Function1.H:92
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
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:346
label size() const noexcept
Definition UList.H:706
void size(const label n)
Definition UList.H:118
const cyclicFvPatch & cyclicPatch() const
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:133
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
virtual void setJump(const Field< vector > &jump)
fixedJumpFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
const fvBoundaryMesh & boundaryMesh() const noexcept
Return boundaryMesh reference.
Definition fvPatch.H:268
const vectorField & Cf() const
Return face centres.
Definition fvPatch.C:113
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
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...
This boundary condition provides a jump condition for U across a cyclic pressure jump condition and a...
virtual tmp< fvPatchField< vector > > clone() const
Return a clone.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
swirlFanVelocityFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal 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:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
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.
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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)
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
fvPatchField< vector > fvPatchVectorField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
Unit conversion functions.