Loading...
Searching...
No Matches
propellerInfo.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) 2021-2025 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
28#include "propellerInfo.H"
29#include "cylindricalCS.H"
30#include "fvMesh.H"
31#include "IOMRFZoneList.H"
33#include "interpolation.H"
34#include "Function1.H"
35#include "surfaceWriter.H"
36#include "treeDataCell.H"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
47}
48}
49
52({
53 { rotationMode::SPECIFIED, "specified" },
54 { rotationMode::MRF, "MRF" },
55});
56
57
58// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
59
61(
62 const dictionary& dict
63)
64{
65 vector origin(Zero);
66 vector axis(Zero);
67
68 switch (rotationMode_)
69 {
71 {
72 origin = dict.get<vector>("origin");
73 axis = dict.get<vector>("axis");
74 axis.normalise();
75
76 // Can specify 'rpm' or 'n' (rev/s)
77 if (dict.readIfPresent("rpm", n_))
78 {
79 n_ /= 60; // -> rev/s
80 dict.readIfPresent("n", n_); // Optional if rpm was specified
81 }
82 else
83 {
84 n_ = dict.get<scalar>("n");
85 }
86 break;
87 }
89 {
90 MRFName_ = dict.get<word>("MRF");
91 const auto* MRFZones =
92 mesh_.cfindObject<IOMRFZoneList>("MRFProperties");
93
94 if (!MRFZones)
95 {
97 << "Unable to find MRFProperties in the database. "
98 << "Is this an MRF case?"
100 }
101 const auto& mrf = MRFZones->MRFZoneList::getFromName(MRFName_);
102
103 dict.readIfPresent("originOffset", origin);
104 origin += mrf.origin();
105 axis = mrf.axis();
106
107 // Convert rad/s to rev/s
108 n_ = (mrf.Omega() & axis)/constant::mathematical::twoPi;
109 break;
110 }
111 default:
112 {
114 << "Unhandled enumeration " << rotationModeNames_[rotationMode_]
115 << abort(FatalError);
116 }
117 }
118
119
120 // Optional orientation axis for cylindrical coordinate system
121 vector alphaAxis;
122 if (dict.readIfPresent("alphaAxis", alphaAxis))
123 {
124 alphaAxis.normalise();
125 coordSysPtr_.reset
126 (
127 new coordSystem::cylindrical(origin, axis, alphaAxis)
128 );
129 }
130 else
132 // Use best-guess for an orthogonal second axis
133 coordSysPtr_.reset(new coordSystem::cylindrical(origin, axis));
134 }
135}
136
137
139{
140 switch (rotationMode_)
141 {
142 case rotationMode::SPECIFIED:
143 {
144 // Set on dictionary re-read
145 break;
146 }
147 case rotationMode::MRF:
148 {
149 const auto* MRFZones =
150 mesh_.cfindObject<IOMRFZoneList>("MRFProperties");
151
152 if (!MRFZones)
153 {
155 << "Unable to find MRFProperties in the database. "
156 << "Is this an MRF case?"
157 << exit(FatalError);
158 }
159
160 const auto& mrf = MRFZones->MRFZoneList::getFromName(MRFName_);
161
162 // Convert rad/s to revolutions per second
163 n_ = (mrf.Omega() & mrf.axis())/constant::mathematical::twoPi;
164 break;
165 }
166 default:
167 {
169 << "Unhandled enumeration " << rotationModeNames_[rotationMode_]
170 << abort(FatalError);
171 }
172 }
173}
174
175
177{
178 if (!writeToFile())
179 {
180 return;
181 }
182
183 if (writePropellerPerformance_ && !propellerPerformanceFilePtr_)
184 {
185 propellerPerformanceFilePtr_ =
186 newFileAtStartTime("propellerPerformance");
187 auto& os = propellerPerformanceFilePtr_();
188
189 writeHeader(os, "Propeller performance");
190 writeHeaderValue(os, "CofR", coordSysPtr_->origin());
191 writeHeaderValue(os, "Radius", radius_);
192 writeHeaderValue(os, "Axis", coordSysPtr_->e3());
193
194 writeHeader(os, "");
195
196 writeCommented(os, "Time");
197 writeTabbed(os, "n");
198 writeTabbed(os, "URef");
199 writeTabbed(os, "J");
200 writeTabbed(os, "KT");
201 writeTabbed(os, "10*KQ");
202 writeTabbed(os, "eta0");
203 os << nl;
204 }
205
206 if (writeWakeFields_)
207 {
210 newFileAtStartTime("axialWake");
211 }
212}
213
214
216{
217 const auto* UPtr = mesh_.cfindObject<volVectorField>(UName_);
218
219 if (!UPtr)
220 {
222 << "Unable to find velocity field " << UName_
223 << " . Available vector fields are: "
224 << flatOutput(mesh_.sortedNames<volVectorField>())
225 << exit(FatalError);
226
228 }
229
230 return *UPtr;
231}
232
233
235(
236 const coordinateSystem& coordSys,
237 const scalar r1,
238 const scalar r2,
239 const scalar nTheta,
240 const label nRadius,
241 faceList& faces,
243) const
244{
245 label nPoint = nRadius*nTheta;
246 if (r1 < SMALL)
247 {
248 nPoint += 1; // 1 for origin
249 }
250 else
251 {
252 nPoint += nTheta;
253 }
254 const label nFace = nRadius*nTheta;
255
256 points.resize_nocopy(nPoint);
257 faces.resize_nocopy(nFace);
258
259 const point& origin = coordSys.origin();
260 const scalar zCoord = 0;
261 label pointi = 0;
262 for (int radiusi = 0; radiusi <= nRadius; ++radiusi)
263 {
264 if (r1 < SMALL && radiusi == 0)
265 {
266 points[pointi++] = origin;
267 }
268 else
269 {
270 const scalar r = r1 + radiusi*(r2 - r1)/nRadius;
271
272 for (label i = 0; i < nTheta; ++i)
273 {
274 point p
275 (
276 r,
277 (i/scalar(nTheta))*constant::mathematical::twoPi,
278 zCoord
279 );
280
281 points[pointi++] = coordSys.globalPosition(p);
282 }
283 }
284 }
285
286
287 const List<label> ptIDs(identity(nTheta));
288
289 // Faces
290 label facei = 0;
291 label pointOffset0 = 0;
292 label radiusOffset = 0;
293 DynamicList<label> facePts(4);
294 for (int radiusi = 0; radiusi < nRadius; ++radiusi)
295 {
296 if (r1 < SMALL && radiusi == 0)
297 {
298 radiusOffset = 1;
299 pointOffset0 = 1;
300
301 // Adding faces as a fan about the origin
302 for (label thetai = 1; thetai <= nTheta; ++thetai)
303 {
304 facePts.clear();
305
306 // Append triangle about origin
307 facePts.append(0);
308 facePts.append(thetai);
309 facePts.append(1 + ptIDs.fcIndex(thetai - 1));
310
311 faces[facei++] = face(facePts);
312 }
313 }
314 else
315 {
316 for (label thetai = 0; thetai < nTheta; ++thetai)
317 {
318 facePts.clear();
319
320 label offset = pointOffset0 + (radiusi-radiusOffset)*nTheta;
321
322 // Inner
323 facePts.append(offset + ptIDs.fcIndex(thetai - 1));
324 facePts.append(offset + ptIDs.fcIndex(thetai));
325
326 // Outer
327 facePts.append(offset + nTheta + ptIDs.fcIndex(thetai));
328 facePts.append(offset + nTheta + ptIDs.fcIndex(thetai - 1));
329
330 faces[facei++] = face(facePts);
331 }
332 }
333 }
334}
335
336
338(
339 const dictionary& dict
340)
341{
342 const dictionary& sampleDiskDict(dict.subDict("sampleDisk"));
343
344 const scalar r1 = sampleDiskDict.getScalar("r1");
345 const scalar r2 = sampleDiskDict.getScalar("r2");
346
347 nTheta_ = sampleDiskDict.getLabel("nTheta");
348 nRadial_ = sampleDiskDict.getLabel("nRadial");
349
350 setSampleDiskGeometry
351 (
352 coordSysPtr_(),
353 r1,
354 r2,
355 nTheta_,
356 nRadial_,
357 faces_,
358 points_
359 );
360
361 // Surface writer (keywords: surfaceWriter, writeOptions)
362
363 word writerType;
364 if (sampleDiskDict.readIfPresent("surfaceWriter", writerType))
365 {
366 surfaceWriterPtr_ = surfaceWriter::New
367 (
368 writerType,
370 (
371 sampleDiskDict,
372 writerType,
373 "writeOptions"
374 )
375 );
376
377 // Use outputDir/TIME/surface-name
378 surfaceWriterPtr_->useTimeDir(true);
379 }
380
382 sampleDiskDict.getOrDefault("errorOnPointNotFound", false);
383
385}
386
387
389{
390 if (!writeWakeFields_)
391 {
392 return;
393 }
394
395 treeBoundBox bb(points_);
396 bb.inflate(0.05);
397 DynamicList<label> treeCellIDs(10*points_.size());
398
399 const auto& meshCells = mesh_.cells();
400 const auto& meshFaces = mesh_.faces();
401 const auto& meshPoints = mesh_.points();
402
403 forAll(meshCells, celli)
404 {
405 bool found = false;
406
407 for (const label facei : meshCells[celli])
408 {
409 for (const label fpi : meshFaces[facei])
410 {
411 if (bb.contains(meshPoints[fpi]))
412 {
413 found = true;
414 break;
415 }
416 }
417
418 if (found)
419 {
420 treeCellIDs.append(celli);
421 break;
422 }
423 }
424 }
425
426 indexedOctree<treeDataCell> tree
427 (
428 treeDataCell(true, mesh_, std::move(treeCellIDs), polyMesh::CELL_TETS),
429 bb,
430 10,
431 100,
432 10
433 );
434
435 cellIds_.setSize(points_.size(), -1);
436 pointMask_.setSize(points_.size(), false);
437
438 // Kick the tet base points calculation to avoid parallel comms later
439 (void)mesh_.tetBasePtIs();
440
441 const auto& treeData = tree.shapes();
442
443 forAll(points_, pointi)
444 {
445 const vector& pos = points_[pointi];
446
447// label meshCelli = mesh_.findCell(pos);
448 label treeCelli = tree.findInside(pos);
449
450 label proci = treeCelli >= 0 ? Pstream::myProcNo() : -1;
451
452 reduce(proci, maxOp<label>());
453
454 pointMask_[pointi] = treeCelli != -1;
455
456 if (proci >= 0)
457 {
458 cellIds_[pointi] =
459 (
460 proci == Pstream::myProcNo()
461 ? treeData.objectIndex(treeCelli)
462 : -1
463 );
464 }
465 else
466 {
467 if (errorOnPointNotFound_)
468 {
470 << "Position " << pos << " not found in mesh"
471 << abort(FatalError);
472 }
473 else
474 {
476 << "Position " << pos << " not found in mesh"
477 << endl;
478 }
480 }
481
483}
484
485
487(
488 const scalarField& field
489) const
490{
491 if (field.size() != points_.size())
492 {
494 << "Inconsistent field sizes: input:" << field.size()
495 << " points:" << points_.size()
496 << abort(FatalError);
497 }
498
499 scalar sumArea = 0;
500 scalar sumFieldArea = 0;
501
502 for (const face& f : faces_)
503 {
504 bool valid = true;
505 scalar faceValue = 0;
506 for (const label pti : f)
507 {
508 // Exclude contributions where sample cell for point was not found
509 if (!pointMask_[pti])
510 {
511 valid = false;
512 break;
513 }
514 faceValue += field[pti];
515 }
516
517 if (valid)
518 {
519 scalar area = f.mag(points_);
520 sumArea += area;
521 sumFieldArea += faceValue/f.size()*area;
523 }
524
525 return sumFieldArea/(sumArea + ROOTVSMALL);
526}
527
528
530(
531 const vectorField& U,
532 const vectorField& Ur,
533 const scalar URef
534)
535{
536 // Write surface
537 if (!surfaceWriterPtr_)
538 {
539 return;
540 }
541
542
543 // Time-aware, with time spliced into the output path
544 surfaceWriterPtr_->isPointData(true);
545 surfaceWriterPtr_->beginTime(time_);
546 surfaceWriterPtr_->open
547 (
548 points_,
549 faces_,
550 (baseFileDir() / name() / "surfaces" / "disk"),
551 false // serial - already merged
552 );
553 surfaceWriterPtr_->nFields(4); // Legacy VTK
554 surfaceWriterPtr_->write("U", U);
555 surfaceWriterPtr_->write("Ur", Ur);
556 surfaceWriterPtr_->write("UNorm", U/URef);
557 surfaceWriterPtr_->write("UrNorm", Ur/URef);
558 surfaceWriterPtr_->endTime();
559 surfaceWriterPtr_->clear();
560}
561
562
564{
565 if (!writePropellerPerformance_)
566 {
567 return;
568 }
569
570 // Update n_
571 setRotationalSpeed();
572
573 const vector sumForce = forceEff();
574 const vector sumMoment = momentEff();
575
576 const scalar diameter = 2*radius_;
577 const scalar URef = URefPtr_->value(time_.timeOutputValue());
578 const scalar j = -URef/mag(n_ + ROOTVSMALL)/diameter;
579 const scalar denom = rhoRef_*sqr(n_)*pow4(diameter);
580 const scalar kt = (sumForce & coordSysPtr_->e3())/denom;
581 const scalar kq =
582 -sign(n_)*(sumMoment & coordSysPtr_->e3())/(denom*diameter);
583 const scalar etaO = kt*j/(kq*constant::mathematical::twoPi + ROOTVSMALL);
584
585 if (writeToFile())
586 {
587 auto& os = propellerPerformanceFilePtr_();
588
589 writeCurrentTime(os);
590 os << tab << n_
591 << tab << URef
592 << tab << j
593 << tab << kt
594 << tab << 10*kq
595 << tab << etaO
596 << nl;
597
598 os.flush();
599 }
600
601 Log << type() << " " << name() << " output:" << nl
602 << " Revolutions per second, n : " << n_ << nl
603 << " Reference velocity, URef : " << URef << nl
604 << " Advance coefficient, J : " << j << nl
605 << " Thrust coefficient, Kt : " << kt << nl
606 << " Torque coefficient, 10*Kq : " << 10*kq << nl
607 << " Efficiency, etaO : " << etaO << nl
608 << nl;
609
610
611 // Write state/results information
612 setResult("n", n_);
613 setResult("URef", URef);
614 setResult("Kt", kt);
615 setResult("Kq", kq);
616 setResult("J", j);
617 setResult("etaO", etaO);
618}
619
620
622(
623 const vectorField& U,
624 const scalar URef
625)
626{
627 if (!Pstream::master()) return;
628
629 // Velocity
630 auto& os = wakeFilePtr_();
631
632 const pointField propPoints(coordSysPtr_->localPosition(points_));
633 const label offset =
634 mag(propPoints[1][0] - propPoints[0][0]) < SMALL ? 0 : 1;
635 const scalar rMax = propPoints.last()[0];
636
637 const scalar UzMean = meanSampleDiskField(U.component(2));
638
639 writeHeaderValue(os, "Time", time_.timeOutputValue());
640 writeHeaderValue(os, "Reference velocity", URef);
641 writeHeaderValue(os, "Direction", coordSysPtr_->e3());
642 writeHeaderValue(os, "Wake", 1 - UzMean/URef);
643 writeHeader(os, "");
644 writeCommented(os, "r/R");
645 writeTabbed(os, "alpha");
646 writeTabbed(os, "(x y z)");
647 writeTabbed(os, "(Ur Utheta Uz)");
648 os << nl;
649
650 for (label thetai = 0; thetai < nTheta_; ++thetai)
651 {
652 const scalar deg = 360*thetai/scalar(nTheta_);
653
654 for (label radiusi = 0; radiusi <= nRadial_; ++radiusi)
655 {
656 label pointi = radiusi*nTheta_ + thetai + offset;
657
658 if (radiusi == 0 && offset == 1)
659 {
660 // Only a single point at the centre - repeat for all thetai
661 pointi = 0;
662 }
663
664 if (pointMask_[pointi])
665 {
666 const scalar rR = propPoints[radiusi*nTheta_][0]/rMax;
667
668 os << rR << tab << deg << tab
669 << points_[pointi] << tab << U[pointi] << nl;
670 }
671 }
672 }
674 writeBreak(os);
675
676 os << endl;
677}
678
679
681(
682 const vectorField& U,
683 const scalar URef
684)
685{
686 if (!Pstream::master()) return;
687
688 // Alternative common format - axial wake component
689 auto& os = axialWakeFilePtr_();
690
691 const pointField propPoints(coordSysPtr_->localPosition(points_));
692 const label offset =
693 mag(propPoints[1][0] - propPoints[0][0]) < SMALL ? 0 : 1;
694 const scalar rMax = propPoints.last()[0];
695
696 writeHeaderValue(os, "Time", time_.timeOutputValue());
697
698 os << "# angle";
699 for (label radiusi = 0; radiusi <= nRadial_; ++radiusi)
700 {
701 label pointi = radiusi*nTheta_;
702 scalar r = propPoints[pointi][0];
703 os << tab << "r/R=" << r/rMax;
704 }
705 os << nl;
706
707 for (label thetai = 0; thetai < nTheta_; ++thetai)
708 {
709 os << 360*thetai/scalar(nTheta_);
710
711 for (label radiusi = 0; radiusi <= nRadial_; ++radiusi)
712 {
713 label pointi = radiusi*nTheta_ + thetai + offset;
714
715 if (radiusi == 0 && offset == 1)
716 {
717 // Only a single point at the centre - repeat for all thetai
718 pointi = 0;
719 }
720
721 if (pointMask_[pointi])
722 {
723 os << tab << 1 - U[pointi][2]/URef;
724 }
725 else
726 {
727 os << tab << "undefined";
728 }
729 }
730
731 os << nl;
732 }
734 writeBreak(os);
735
736 os << endl;
737}
738
739
741{
742 if (!writeWakeFields_)
743 {
744 return;
745 }
746
747 scalar URef = URef0;
748 if (mag(URef) < ROOTSMALL)
749 {
751 << "Magnitude of reference velocity should be greater than zero"
752 << endl;
753
754 URef = ROOTVSMALL;
755 }
756
757 // Normalised velocity
758 const vectorField UDisk(interpolate(U(), vector::uniform(nanValue_))());
759 const vectorField UrDisk(coordSysPtr_->localVector(UDisk));
760
761 // Surface field data
762 writeSampleDiskSurface(UDisk, UrDisk, URef);
763
764 // Write wake text files
765 writeWake(UrDisk, URef);
766 writeAxialWake(UrDisk, URef);
767}
768
769
770template<class Type>
772(
774 const Type& defaultValue
775) const
776{
777 auto tfield = tmp<Field<Type>>::New(points_.size(), defaultValue);
778 auto& field = tfield.ref();
779
780 autoPtr<interpolation<Type>> interpolator
781 (
782 interpolation<Type>::New(interpolationScheme_, psi)
783 );
784
785 forAll(points_, pointi)
786 {
787 const label celli = cellIds_[pointi];
788
789 if (cellIds_[pointi] != -1)
790 {
791 const point& position = points_[pointi];
792 field[pointi] = interpolator().interpolate(position, celli);
793 }
794 }
795
796 Pstream::listReduce(field, maxOp<Type>());
798 return tfield;
799}
800
801
802// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
803
805(
806 const word& name,
807 const Time& runTime,
808 const dictionary& dict,
809 bool readFields
810)
811:
812 forces(name, runTime, dict, false),
813 dict_(dict),
814 radius_(0),
815 URefPtr_(nullptr),
816 rotationMode_(rotationMode::SPECIFIED),
817 MRFName_(),
818 n_(0),
819 writePropellerPerformance_(true),
820 propellerPerformanceFilePtr_(nullptr),
821 writeWakeFields_(true),
822 surfaceWriterPtr_(nullptr),
823 nTheta_(0),
824 nRadial_(0),
825 points_(),
826 errorOnPointNotFound_(false),
827 faces_(),
828 cellIds_(),
829 pointMask_(),
830 interpolationScheme_("cell"),
831 wakeFilePtr_(nullptr),
832 axialWakeFilePtr_(nullptr),
833 nanValue_(pTraits<scalar>::min),
834 initialised_(false)
835{
836 if (readFields)
838 read(dict);
839 Log << endl;
840 }
841}
842
843
845(
846 const word& name,
847 const objectRegistry& obr,
848 const dictionary& dict,
849 bool readFields
850)
852 propellerInfo(name, obr.time(), dict, false)
853{}
854
855
856// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
857
859{
860 if (forces::read(dict))
861 {
862 dict_ = dict;
863
864 radius_ = dict.getScalar("radius");
865 URefPtr_.reset(Function1<scalar>::New("URef", dict, &mesh_));
866 rotationMode_ = rotationModeNames_.get("rotationMode", dict);
867
868 writePropellerPerformance_ =
869 dict.get<bool>("writePropellerPerformance");
870
871 writeWakeFields_ = dict.get<bool>("writeWakeFields");
872 if (writeWakeFields_)
873 {
874 dict.readIfPresent("interpolationScheme", interpolationScheme_);
875
876 dict.readIfPresent("nanValue", nanValue_);
877 }
878
879 return true;
880 }
881
882 return false;
883}
884
885
887{
888 if (!initialised_)
889 {
890 // Must be set before setting the surface
891 setCoordinateSystem(dict_);
892
893 if (writeWakeFields_)
894 {
895 setSampleDiskSurface(dict_);
896 }
897
898 initialised_ = true;
899 }
900
901 calcForcesMoments();
902
903 createFiles();
904
905 if (writeWakeFields_)
906 {
907 // Only setting mean axial velocity result during execute
908 // - wake fields are 'heavy' and controlled separately using the
909 // writeControl
910 const vectorField
911 UDisk
912 (
913 coordSysPtr_->localVector
914 (
916 (
917 U(),
918 vector::uniform(nanValue_)
919 )()
920 )
921 );
922 const scalar UzMean = meanSampleDiskField(UDisk.component(2));
923
924 setResult("UzMean", UzMean);
925 }
928
929 return true;
930}
931
932
934{
935 const scalar URef = URefPtr_->value(time_.timeOutputValue());
936 writeWakeFields(URef);
937
938 return true;
939}
940
945}
946
947
949{
950 updateSampleDiskCells();
951}
952
953
954// ************************************************************************* //
#define Log
Definition PDRblock.C:28
bool found
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
void append(const T &val)
Copy append an element to the end of this list.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
Generic GeometricField class.
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries,...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
static void listCombineReduce(UList< T > &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::listReduce with an in-place cop.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & last()
Access last element of the list, position [size()-1].
Definition UList.H:971
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition VectorI.H:114
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition boundBoxI.H:381
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Base class for coordinate system specification, the default coordinate system type is cartesian .
point globalPosition(const point &local) const
From local coordinate position to global (cartesian) position.
virtual const point & origin() const
Return origin.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
label getLabel(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<label>(const word&, keyType::option).
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<scalar>(const word&, keyType::option).
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Abstract base-class for Time/database function objects.
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition forces.H:321
virtual void calcForcesMoments()
Calculate forces and moments.
Definition forces.C:665
scalar rhoRef_
Reference density needed for incompressible calculations.
Definition forces.H:387
word UName_
Name of velocity field.
Definition forces.H:402
forces(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
Definition forces.C:510
virtual vector forceEff() const
Return the total force.
Definition forces.C:783
autoPtr< coordinateSystem > coordSysPtr_
Coordinate system used when evaluating forces and moments.
Definition forces.H:377
virtual vector momentEff() const
Return the total moment.
Definition forces.C:789
virtual bool read(const dictionary &)
Read the function-object dictionary.
Definition forces.C:589
const fvMesh & mesh_
Reference to the fvMesh.
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).
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.
dictionary dict_
Copy of dictionary used during construction.
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.
bool errorOnPointNotFound_
Flag to raise an error if the sample point is not found in the mesh. Default = false to enable....
void setCoordinateSystem(const dictionary &dict)
Set the coordinate system.
scalar meanSampleDiskField(const scalarField &field) const
Return the area average of a field.
void setRotationalSpeed()
Set the rotational speed.
pointField points_
Surface points.
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.
scalar nanValue_
Default value when a sample point is not found; default = scalar::min.
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 the function-object operations.
virtual bool write()
Write the function-object results.
void updateSampleDiskCells()
Set the sample cells corresponding to the sample points.
bool writeWakeFields_
Flag to write wake fields.
static const Enum< rotationMode > rotationModeNames_
virtual bool read(const dictionary &)
Read the function-object dictionary.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
void setResult(const word &entryName, const Type &value)
Add result.
const Time & time_
Reference to the time database.
const Time & time() const
Return time database.
fileName baseFileDir() const
Return the base directory for output.
Definition writeFile.C:43
virtual autoPtr< OFstream > newFileAtStartTime(const word &name) const
Return autoPtr to a new file using the simulation start time.
Definition writeFile.C:156
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
virtual void writeBreak(Ostream &os) const
Write a break marker to the stream.
Definition writeFile.C:367
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
virtual bool writeToFile() const
Flag to allow writing to file.
Definition writeFile.C:286
Non-pointer based hierarchical recursive searching.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
A class for managing temporary objects.
Definition tmp.H:75
Standard boundBox with extra functionality for use in octree.
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
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
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
const volScalarField & psi
rDeltaTY field()
dynamicFvMesh & mesh
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar twoPi(2 *M_PI)
const wordList area
Standard area field types (scalar, vector, tensor, etc).
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static void writeHeader(Ostream &os, const word &fieldName)
List< face > faceList
List of faces.
Definition faceListFwd.H:41
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
dimensionedScalar pow4(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
static bool initialised_(false)
labelList f(nPoints)
dictionary dict
Tree tree(triangles.begin(), triangles.end())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299