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