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-2025 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace fv
34 {
37  (
38  fvModel,
41  );
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 {
50  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
51 
52  UName_ = dict.lookupOrDefault<word>
53  (
54  "U",
56  );
57 
58  if (dict.found("centre"))
59  {
60  centre_ = dict.lookup<vector>("centre");
61  }
62  else
63  {
64  const Field<vector> zoneCellCentres(mesh().cellCentres(), zone_.zone());
65  const Field<scalar> zoneCellVolumes(mesh().cellVolumes(), zone_.zone());
66  centre_ = gSum(zoneCellVolumes*zoneCellCentres)/zone_.V();
67  }
68 
69  normal_ = normalised(dict.lookup<vector>("normal"));
70 
71  n_ = dict.lookup<scalar>("n");
72  rotationDir_ = sign(n_);
73  n_ = mag(n_);
74 
75  dProp_ = dict.lookup<scalar>("dPropeller");
76  dHub_ = dict.lookup<scalar>("dHub");
77 
79  Function1<vector2D>::New("propellerCurve", dimless, dimless, dict);
80 
81  log_ = dict.lookupOrDefault<Switch>("log", false);
82 
83  if (log_ && Pstream::master())
84  {
85  logFile_ = new functionObjects::logFile(mesh(), name(), typeName);
86  logFile_->writeTimeColumnHeaders
87  (
88  {
89  "n",
90  "J",
91  "Jcorr",
92  "Udisk",
93  "Ucorr",
94  "Kt",
95  "Kq",
96  "T/rho",
97  "Q/rho",
98  "force",
99  "moment"
100  }
101  );
102  }
103 }
104 
105 
106 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
107 
108 Foam::scalar Foam::fv::propellerDisk::diskThickness(const vector& centre) const
109 {
110  const Field<vector> zoneCellCentres(mesh().cellCentres(), zone_.zone());
111  const scalar r2 = gMax(magSqr(zoneCellCentres - centre));
112 
113  return zone_.V()/(constant::mathematical::pi*r2);
114 }
115 
116 
118 (
119  const vectorField& U,
120  const vector& nHat
121 ) const
122 {
123  const labelList& cells = zone_.zone();
124  const scalarField& V = mesh().V();
125 
126  scalar VUn = 0;
127  forAll(cells, i)
128  {
129  VUn += V[cells[i]]*(nHat & U[cells[i]]);
130  }
131  reduce(VUn, sumOp<scalar>());
132 
133  const scalar Uref = VUn/zone_.V();
134 
135  return mag(Uref/(n_*dProp_));
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
140 
142 (
143  const word& name,
144  const word& modelType,
145  const fvMesh& mesh,
146  const dictionary& dict
147 )
148 :
149  fvModel(name, modelType, mesh, dict),
150  zone_(mesh, dict),
151  phaseName_(word::null),
152  UName_(word::null),
153  force_(Zero),
154  moment_(Zero),
155  log_(false)
156 {
158 }
159 
160 
161 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162 
164 {
165  if (fvModel::read(dict))
166  {
167  zone_.read(coeffs(dict));
168  readCoeffs(coeffs(dict));
169  return true;
170  }
171  else
172  {
173  return false;
174  }
175 }
176 
177 
179 {
180  return wordList(1, UName_);
181 }
182 
183 
185 (
186  const volScalarField& alpha,
187  const volScalarField& rho,
188  const volVectorField& U,
189  fvMatrix<vector>& eqn
190 ) const
191 {
192  addActuationDiskAxialInertialResistance
193  (
194  eqn.source(),
195  alpha,
196  rho,
197  U
198  );
199 }
200 
201 
203 (
204  const volVectorField& U,
205  fvMatrix<vector>& eqn
206 ) const
207 {
208  addActuationDiskAxialInertialResistance
209  (
210  eqn.source(),
213  U
214  );
215 }
216 
217 
219 (
220  const volScalarField& rho,
221  const volVectorField& U,
222  fvMatrix<vector>& eqn
223 ) const
224 {
225  addActuationDiskAxialInertialResistance
226  (
227  eqn.source(),
229  rho,
230  U
231  );
232 }
233 
234 
236 {
237  zone_.movePoints();
238  return true;
239 }
240 
241 
243 {
244  zone_.topoChange(map);
245 }
246 
247 
249 {
250  zone_.mapMesh(map);
251 }
252 
253 
255 (
256  const polyDistributionMap& map
257 )
258 {
259  zone_.distribute(map);
260 }
261 
262 
263 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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)
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
functionObject support for writing log files
Definition: logFile.H:60
scalar V() const
Return const access to the total cell volume.
Definition: fvCellZoneI.H:34
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:96
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
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.
vector centre_
Propeller disk centre.
word UName_
Name of the velocity field.
virtual void readCoeffs(const dictionary &dict)
Read the model coefficients.
Definition: propellerDisk.C:48
word phaseName_
The name of the phase to which this fvModel applies.
fvCellZone zone_
The cellZone the fvConstraint applies to.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
vector normal_
Propeller disk normal direction.
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.
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.
const cellZone & zone() const
Return const access to the cell set.
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Type gMax(const FieldField< Field, Type > &f)
labelList fv(nPoints)
dictionary dict