propellerInfo.H
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) 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
26Class
27 Foam::functionObjects::propellerInfo
28
29Group
30 grpForcesFunctionObjects
31
32Description
33 Calculates propeller performance and wake field properties.
34
35 Controlled by executeControl:
36 - Propeller performance
37 - Thrust coefficient, Kt
38 - Torque coefficient, 10*Kq
39 - Advance coefficient, J
40 - Open water efficiency, etaO
41 - Written to postProcessing/<name>/<time>/propellerPerformance.dat
42
43 Controlled by writeControl:
44 - Wake field text file
45 - Wake: 1 - UzMean/URef
46 - Velocity in cylindrical coordinates at xyz locations
47 - Written to postProcessing/<name>/<time>/wake.dat
48 - Axial wake field text file
49 - 1 - Uz/URef at r/R and angle
50 - Written to postProcessing/<name>/<time>/axialWake.dat
51 - Velocity surface
52 - Written to postProcessing/<name>/surfaces/<time>/disk.<fileType>
53
54Usage
55 Example of function object specification:
56 \verbatim
57 propellerInfo1
58 {
59 type propellerInfo;
60 libs (forces);
61 writeControl writeTime;
62
63 patches ("propeller.*");
64
65 URef 5; // Function1 type; 'constant' form shown here
66 rho rhoInf; // incompressible
67 rhoInf 1.2;
68
69 // Optionally write propeller performance data
70 writePropellerPerformance yes;
71
72
73 // Propeller data:
74
75 // Radius
76 radius 0.1;
77
78 rotationMode specified; // specified | MRF
79
80 // rotationMode = specified:
81 origin (0 -0.1 0);
82 n 25.15;
83 axis (0 1 0);
84
85 // Optional reference direction for angle (alpha) = 0
86 alphaAxis (1 0 0);
87
88
92
93 // Optionally write wake text files
94 // Note: controlled by writeControl
95 writeWakeFields yes;
96
97 // Sample plane (disk) properties
98 // Note: controlled by writeControl
99 sampleDisk
100 {
101 surfaceWriter vtk;
102 r1 0.05;
103 r2 0.2;
104 nTheta 36;
105 nRadial 10;
106 interpolationScheme cellPoint;
107 errorOnPointNotFound false;
108 }
109 }
110 \endverbatim
111
112 Where the entries comprise:
113 \table
114 Property | Description | Required | Deflt value
115 type | Type name: propellerInfo | yes |
116 log | Write to standard output | no | no
117 patches | Patches included in the forces calculation | yes |
118 p | Pressure field name | no | p
119 U | Velocity field name | no | U
120 rho | Density field name | no | rho
121 URef | Reference velocity | yes |
122 rotationMode | Rotation mode (see below) | yes |
123 origin | Sample disk centre | no* |
124 n | Revolutions per second | no* |
125 axis | Propeller axis | no* |
126 alphaAxis | Axis that defines alpha=0 dir | no |
127 MRF | Name of MRF zone | no* |
128 originOffset | Origin offset for MRF mode | no | (0 0 0)
129 writePropellerPerformance| Write propeller performance text file | yes |
130 writeWakeFields | Write wake field text files | yes |
131 surfaceWriter | Sample disk surface writer | no* |
132 r1 | Sample disk inner radius | no | 0
133 r2 | Sample disk outer radius | no* |
134 nTheta | Divisions in theta direction | no* |
135 nRadial | Divisions in radial direction | no* |
136 interpolationScheme | Sampling interpolation scheme | no* | cell
137 \endtable
138
139
140Note
141- URef is a scalar Function1 type, i.e. supports constant, table, lookup values
142- rotationMode is used to set the origin, axis and revolutions per second
143 - if set to 'specified' all 3 entries are required
144 - note: origin is the sample disk origin
145 - if set to 'MRF' only the MRF entry is required
146 - to move the sample disk away from the MRF origin, use the originOffset
147- if writePropellerPerformance is set to on|true:
148 - propellerPerformance text file will be written
149- if writeWakeFields is set to on|true:
150 - wake and axialWake text files will be written
151- if the surfaceWriter entry is set, the sample disk surface will be written
152 - extents set according to the r1 and r2 entries
153 - discretised according to the nTheta and nRadial entries
154
155See also
156 Foam::functionObject::forces
157 Foam::Function1
158
159SourceFiles
160 propellerInfo.C
161
162\*---------------------------------------------------------------------------*/
163
164#ifndef functionObjects_propellerInfo_H
165#define functionObjects_propellerInfo_H
166
167#include "forces.H"
168#include "Enum.H"
169#include "faceList.H"
170
171// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172
173namespace Foam
174{
175
176template<class Type>
177class Function1;
178
179class surfaceWriter;
180
181namespace functionObjects
182{
183
184/*---------------------------------------------------------------------------*\
185 Class propellerInfo Declaration
186\*---------------------------------------------------------------------------*/
187
188class propellerInfo
189:
190 public forces
191{
192
193public:
194
195 enum class rotationMode
196 {
197 SPECIFIED,
198 MRF
199 };
200
201 static const Enum<rotationMode> rotationModeNames_;
202
203
204protected:
205
206 // Protected data
207
208 //- Propeller radius
209 scalar radius_;
210
211 //- Reference velocity
212 autoPtr<Function1<scalar>> URefPtr_;
213
214 //- Rotation mode
216
217 //- Name of MRF zone (if applicable)
218 word MRFName_;
219
220 //- Propeller speed (revolutions per second)
221 scalar n_;
222
223 //- Flag to write performance data
225
226 //- Propeller performance file
227 autoPtr<OFstream> propellerPerformanceFilePtr_;
228
229 //- Flag to write wake fields
230 bool writeWakeFields_;
231
232
233 // Wake field output
234
235 //- Surface writer
236 autoPtr<surfaceWriter> surfaceWriterPtr_;
237
238 //- Number of surface divisions in theta direction
239 label nTheta_;
240
241 //- Number of surface divisions in radial direction
242 label nRadial_;
243
244 //- Surface points
246
247 //- Flag to raise an error if the sample point is not found in the
248 //- mesh. Default = false to enable. e.g. reduced geometry/symmetric
249 //- cases
251
252 //- Surface faces
254
255 //- Surface point sample cell IDs
257
258 //- List of participating points (parallel reduced)
260
261 //- Interpolation scheme
263
264 //- Wake field file
265 autoPtr<OFstream> wakeFilePtr_;
266
267 //- Axial wake field file
268 autoPtr<OFstream> axialWakeFilePtr_;
269
270 //- Default value when a sample point is not found; default =
271 //- scalar::min
272 scalar nanValue_;
273
274
275 // Protected Member Functions
276
277 //- Set the coordinate system
278 void setCoordinateSystem(const dictionary& dict);
279
280 //- Set the rotational speed
281 void setRotationalSpeed();
282
283 //- Create output files
284 void createFiles();
285
286 //- Return the velocity field
287 const volVectorField& U() const;
288
289
290 // Propeller performance text file
291
292 //- Write the wake fields
294
295
296 // Wake text files
297
298 //- Write the wake text file
299 void writeWake(const vectorField& U, const scalar URef);
300
301 //- Write the axial wake text file
302 void writeAxialWake(const vectorField& U, const scalar URef);
303
304 //- Write the wake fields
305 void writeWakeFields(const scalar URef);
306
307
308 // Sample surface functions
310 //- Set the faces and points for the sample surface
312 (
313 const coordinateSystem& coordSys,
314 const scalar r1,
315 const scalar r2,
316 const scalar nTheta,
317 const label nRadius,
318 faceList& faces,
320 ) const;
321
322 //- Set the sample surface based on dictionary settings
324
325 //- Set the sample cells corresponding to the sample points
327
328 //- Return the area average of a field
330
331 //- Write the sample surface
333 (
334 const vectorField& U,
335 const vectorField& Ur,
336 const scalar URef
337 );
339 //- Interpolate from the mesh onto the sample surface
340 template<class Type>
342 (
344 const Type& defaultValue
345 ) const;
346
347 //- No copy construct
348 propellerInfo(const propellerInfo&) = delete;
349
350 //- No copy assignment
351 void operator=(const propellerInfo&) = delete;
352
354public:
355
356 //- Runtime type information
357 TypeName("propellerInfo");
358
360 // Constructors
361
362 //- Construct from Time and dictionary
365 const word& name,
366 const Time& runTime,
368 const bool readFields = true
369 );
371 //- Construct from objectRegistry and dictionary
374 const word& name,
375 const objectRegistry& obr,
377 const bool readFields = true
378 );
380
381 //- Destructor
382 virtual ~propellerInfo() = default;
383
384
385 // Member Functions
387 //- Read the forces data
388 virtual bool read(const dictionary&);
389
390 //- Execute, currently does nothing
391 virtual bool execute();
392
393 //- Write the forces
394 virtual bool write();
395
396 void UpdateMesh(const mapPolyMesh& mpm);
397
398 void movePoints(const polyMesh& mesh);
399};
400
401
402// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
403
404} // End namespace functionObjects
405} // End namespace Foam
406
407// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
408
409#endif
410
411// ************************************************************************* //
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Base class for coordinate system specification, the default coordinate system type is cartesian .
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const word & name() const noexcept
Return the name of this functionObject.
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition: forces.H:325
Calculates propeller performance and wake field properties.
word MRFName_
Name of MRF zone (if applicable)
label nTheta_
Number of surface divisions in theta direction.
boolList pointMask_
List of participating points (parallel reduced)
virtual ~propellerInfo()=default
Destructor.
autoPtr< OFstream > propellerPerformanceFilePtr_
Propeller performance file.
void UpdateMesh(const mapPolyMesh &mpm)
labelList cellIds_
Surface point sample cell IDs.
void writeWake(const vectorField &U, const scalar URef)
Write the wake text file.
void writeWakeFields(const scalar URef)
Write the wake fields.
void createFiles()
Create output files.
rotationMode rotationMode_
Rotation mode.
void setSampleDiskSurface(const dictionary &dict)
Set the sample surface based on dictionary settings.
void writeAxialWake(const vectorField &U, const scalar URef)
Write the axial wake text file.
autoPtr< OFstream > wakeFilePtr_
Wake field file.
void writeSampleDiskSurface(const vectorField &U, const vectorField &Ur, const scalar URef)
Write the sample surface.
const volVectorField & U() const
Return the velocity field.
propellerInfo(const propellerInfo &)=delete
No copy construct.
autoPtr< OFstream > axialWakeFilePtr_
Axial wake field file.
void writePropellerPerformance()
Write the wake fields.
void movePoints(const polyMesh &mesh)
Update for changes of mesh.
tmp< Field< Type > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &psi, const Type &defaultValue) const
Interpolate from the mesh onto the sample surface.
autoPtr< Function1< scalar > > URefPtr_
Reference velocity.
bool writePropellerPerformance_
Flag to write performance data.
void setCoordinateSystem(const dictionary &dict)
Set the coordinate system.
Definition: propellerInfo.C:61
scalar meanSampleDiskField(const scalarField &field) const
Return the area average of a field.
void setRotationalSpeed()
Set the rotational speed.
pointField points_
Surface points.
scalar radius_
Propeller radius.
void setSampleDiskGeometry(const coordinateSystem &coordSys, const scalar r1, const scalar r2, const scalar nTheta, const label nRadius, faceList &faces, pointField &points) const
Set the faces and points for the sample surface.
word interpolationScheme_
Interpolation scheme.
label nRadial_
Number of surface divisions in radial direction.
autoPtr< surfaceWriter > surfaceWriterPtr_
Surface writer.
scalar n_
Propeller speed (revolutions per second)
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Write the forces.
void updateSampleDiskCells()
Set the sample cells corresponding to the sample points.
void operator=(const propellerInfo &)=delete
No copy assignment.
TypeName("propellerInfo")
Runtime type information.
bool writeWakeFields_
Flag to write wake fields.
static const Enum< rotationMode > rotationModeNames_
virtual bool read(const dictionary &)
Read the forces data.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
rDeltaTY field()
dynamicFvMesh & mesh
engineTime & runTime
const pointField & points
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
List< label > labelList
A List of labels.
Definition: List.H:66
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition: List.H:64
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73