propellerDisk.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2024 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 
26 #include "propellerDisk.H"
27 #include "makeFunction1s.H"
28 #include "makeTableReaders.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 namespace fv
37 {
40  (
41  fvModel,
44  );
45 }
46 
49 
50 }
51 
52 
53 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::fv::propellerDisk::readCoeffs(const dictionary& dict)
56 {
57  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
58 
59  UName_ = dict.lookupOrDefault<word>
60  (
61  "U",
63  );
64 
65  diskNormal_ = dict.lookup<vector>("diskNormal");
66  if (mag(diskNormal_) < small)
67  {
69  << "disk direction vector is approximately zero"
70  << exit(FatalIOError);
71  }
72 
73  n_ = dict.lookup<scalar>("n");
74  rotationDir_ = sign(n_);
75  n_ = mag(n_);
76 
77  dProp_ = dict.lookup<scalar>("dPropeller");
78  dHub_ = dict.lookup<scalar>("dHub");
79 
81  Function1<vector2D>::New("propellerCurve", dimless, dimless, dict);
82 
83  log_ = dict.lookupOrDefault<Switch>("log", false);
84 
85  if (log_ && Pstream::master())
86  {
87  logFile_ = new functionObjects::logFile(mesh(), name(), typeName);
88  logFile_->writeTimeColumnHeaders
89  (
90  {"n", "J", "Jcorr", "Udisk", "Ucorr", "Kt", "Kq", "T/rho", "Q/rho"}
91  );
92  }
93 }
94 
95 
96 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
97 
99 {
100  const Field<vector> zoneCellCentres(mesh().cellCentres(), set_.cells());
101  const Field<scalar> zoneCellVolumes(mesh().cellVolumes(), set_.cells());
102 
103  return gSum(zoneCellVolumes*zoneCellCentres)/set_.V();
104 }
105 
106 
107 Foam::scalar Foam::fv::propellerDisk::diskThickness(const vector& centre) const
108 {
109  const Field<vector> zoneCellCentres(mesh().cellCentres(), set_.cells());
110  const scalar r2 = gMax(magSqr(zoneCellCentres - centre));
111 
112  return set_.V()/(constant::mathematical::pi*r2);
113 }
114 
115 
117 (
118  const vectorField& U,
119  const vector& nHat
120 ) const
121 {
122  const labelUList& cells = set_.cells();
123  const scalarField& V = mesh().V();
124 
125  scalar VUn = 0;
126  forAll(cells, i)
127  {
128  VUn += V[cells[i]]*(nHat & U[cells[i]]);
129  }
130  reduce(VUn, sumOp<scalar>());
131 
132  const scalar Uref = VUn/set_.V();
133 
134  return mag(Uref/(n_*dProp_));
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139 
141 (
142  const word& name,
143  const word& modelType,
144  const fvMesh& mesh,
145  const dictionary& dict
146 )
147 :
148  fvModel(name, modelType, mesh, dict),
149  set_(mesh, coeffs()),
150  phaseName_(word::null),
151  UName_(word::null),
152  log_(false)
153 {
154  readCoeffs(coeffs());
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
161 {
162  if (fvModel::read(dict))
163  {
164  set_.read(coeffs());
165  readCoeffs(coeffs());
166  return true;
167  }
168  else
169  {
170  return false;
171  }
172 }
173 
174 
176 {
177  return wordList(1, UName_);
178 }
179 
180 
182 (
183  const volScalarField& alpha,
184  const volScalarField& rho,
185  const volVectorField& U,
186  fvMatrix<vector>& eqn
187 ) const
188 {
189  addActuationDiskAxialInertialResistance
190  (
191  eqn.source(),
192  alpha,
193  rho,
194  U
195  );
196 }
197 
198 
200 (
201  const volVectorField& U,
202  fvMatrix<vector>& eqn
203 ) const
204 {
205  addActuationDiskAxialInertialResistance
206  (
207  eqn.source(),
210  U
211  );
212 }
213 
214 
216 (
217  const volScalarField& rho,
218  const volVectorField& U,
219  fvMatrix<vector>& eqn
220 ) const
221 {
222  addActuationDiskAxialInertialResistance
223  (
224  eqn.source(),
226  rho,
227  U
228  );
229 }
230 
231 
233 {
234  set_.movePoints();
235  return true;
236 }
237 
238 
240 {
241  set_.topoChange(map);
242 }
243 
244 
246 {
247  set_.mapMesh(map);
248 }
249 
250 
252 (
253  const polyDistributionMap& map
254 )
255 {
256  set_.distribute(map);
257 }
258 
259 
260 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
static word groupName(Name name, const word &group)
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:47
Disk momentum source which approximates a propeller based on a given propeller curve.
scalar dProp_
Propeller diameter.
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
word UName_
Name of the velocity field.
word phaseName_
The name of the phase to which this fvModel applies.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
scalar dHub_
Hub diameter.
virtual bool read(const dictionary &dict)
Read dictionary.
scalar diskThickness(const vector &centre) const
Computes the thickness of the disk in streamwise direction.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
autoPtr< Function1< vector2D > > propellerFunction_
Propeller function.
vector diskNormal_
Disk normal direction.
scalar rotationDir_
Rotation direction (obtained from the sign of n_)
scalar J(const vectorField &U, const vector &nHat) const
Return the normalised flow-rate through the disk.
propellerDisk(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
autoPtr< functionObjects::logFile > logFile_
Optional log file.
scalar n_
Rotation speed [1/s].
Switch log_
Optional switch to enable logging of integral properties.
vector diskCentre() const
Computes the disk centre.
Definition: propellerDisk.C:98
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar sign(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
IOerror FatalIOError
makeFoamTableReaders(avType, nullArg)
makeFunction1s(avType, nullArg)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)
labelList fv(nPoints)
dictionary dict