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 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  {
70  case rotationMode::SPECIFIED:
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  }
130 
131  alphaAxis.normalise();
132 
133  coordSysPtr_.reset(new coordSystem::cylindrical(origin, axis, alphaAxis));
134 }
135 
136 
138 {
139  switch (rotationMode_)
140  {
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_ = createFile("propellerPerformance");
185  auto& os = propellerPerformanceFilePtr_();
186 
187  writeHeader(os, "Propeller performance");
188  writeHeaderValue(os, "CofR", coordSysPtr_->origin());
189  writeHeaderValue(os, "Radius", radius_);
190  writeHeaderValue(os, "Axis", coordSysPtr_->e3());
191 
192  writeHeader(os, "");
193 
194  writeCommented(os, "Time");
195  writeTabbed(os, "n");
196  writeTabbed(os, "URef");
197  writeTabbed(os, "J");
198  writeTabbed(os, "KT");
199  writeTabbed(os, "10*KQ");
200  writeTabbed(os, "eta0");
201  os << nl;
202  }
203 
204  if (writeWakeFields_)
205  {
206  if (!wakeFilePtr_) wakeFilePtr_ = createFile("wake");
207  if (!axialWakeFilePtr_) axialWakeFilePtr_ = createFile("axialWake");
208  }
209 }
210 
211 
213 {
214  const auto* UPtr = mesh_.cfindObject<volVectorField>(UName_);
215 
216  if (!UPtr)
217  {
219  << "Unable to find velocity field " << UName_
220  << " . Available vector fields are: "
221  << mesh_.names<volVectorField>()
222  << exit(FatalError);
223  }
224 
225  return *UPtr;
226 }
227 
228 
230 (
231  const coordinateSystem& coordSys,
232  const scalar r1,
233  const scalar r2,
234  const scalar nTheta,
235  const label nRadius,
236  faceList& faces,
238 ) const
239 {
240  label nPoint = nRadius*nTheta;
241  if (r1 < SMALL)
242  {
243  nPoint += 1; // 1 for origin
244  }
245  else
246  {
247  nPoint += nTheta;
248  }
249  const label nFace = nRadius*nTheta;
250 
251  points.setSize(nPoint);
252  faces.setSize(nFace);
253 
254  const point& origin = coordSys.origin();
255  const scalar zCoord = 0;
256  label pointi = 0;
257  for (int radiusi = 0; radiusi <= nRadius; ++radiusi)
258  {
259  if (r1 < SMALL && radiusi == 0)
260  {
261  points[pointi++] = origin;
262  }
263  else
264  {
265  const scalar r = r1 + radiusi*(r2 - r1)/nRadius;
266 
267  for (label i = 0; i < nTheta; ++i)
268  {
269  point p
270  (
271  r,
272  (i/scalar(nTheta))*constant::mathematical::twoPi,
273  zCoord
274  );
275 
276  points[pointi++] = coordSys.globalPosition(p);
277  }
278  }
279  }
280 
281 
282  const List<label> ptIDs(identity(nTheta));
283 
284  // Faces
285  label facei = 0;
286  label pointOffset0 = 0;
287  label radiusOffset = 0;
288  DynamicList<label> facePts(4);
289  for (int radiusi = 0; radiusi < nRadius; ++radiusi)
290  {
291  if (r1 < SMALL && radiusi == 0)
292  {
293  radiusOffset = 1;
294  pointOffset0 = 1;
295 
296  // Adding faces as a fan about the origin
297  for (label thetai = 1; thetai <= nTheta; ++thetai)
298  {
299  facePts.clear();
300 
301  // Append triangle about origin
302  facePts.append(0);
303  facePts.append(thetai);
304  facePts.append(1 + ptIDs.fcIndex(thetai - 1));
305 
306  faces[facei++] = face(facePts);
307  }
308  }
309  else
310  {
311  for (label thetai = 0; thetai < nTheta; ++thetai)
312  {
313  facePts.clear();
314 
315  label offset = pointOffset0 + (radiusi-radiusOffset)*nTheta;
316 
317  // Inner
318  facePts.append(offset + ptIDs.fcIndex(thetai - 1));
319  facePts.append(offset + ptIDs.fcIndex(thetai));
320 
321  // Outer
322  facePts.append(offset + nTheta + ptIDs.fcIndex(thetai));
323  facePts.append(offset + nTheta + ptIDs.fcIndex(thetai - 1));
324 
325  faces[facei++] = face(facePts);
326  }
327  }
328  }
329 }
330 
331 
333 (
334  const dictionary& dict
335 )
336 {
337  const dictionary& sampleDiskDict(dict.subDict("sampleDisk"));
338 
339  const scalar r1 = sampleDiskDict.getScalar("r1");
340  const scalar r2 = sampleDiskDict.getScalar("r2");
341 
342  nTheta_ = sampleDiskDict.getLabel("nTheta");
343  nRadial_ = sampleDiskDict.getLabel("nRadial");
344 
345  setSampleDiskGeometry
346  (
347  coordSysPtr_(),
348  r1,
349  r2,
350  nTheta_,
351  nRadial_,
352  faces_,
353  points_
354  );
355 
356  // Surface writer
357  word surfWriterType;
358  if (sampleDiskDict.readIfPresent("surfaceWriter", surfWriterType))
359  {
360  const auto writeOptions = sampleDiskDict.subOrEmptyDict("writeOptions");
361 
362  surfaceWriterPtr_ = surfaceWriter::New
363  (
364  surfWriterType,
365  writeOptions.subOrEmptyDict(surfWriterType)
366  );
367 
368  // Use outputDir/TIME/surface-name
369  surfaceWriterPtr_->useTimeDir(true);
370  }
371 
372  errorOnPointNotFound_ =
373  sampleDiskDict.getOrDefault("errorOnPointNotFound", false);
374 
375  updateSampleDiskCells();
376 }
377 
378 
380 {
381  if (!writeWakeFields_)
382  {
383  return;
384  }
385 
386  treeBoundBox bb(points_);
387  bb.inflate(0.05);
388  DynamicList<label> treeCellIDs(10*points_.size());
389 
390  const auto& meshCells = mesh_.cells();
391  const auto& meshFaces = mesh_.faces();
392  const auto& meshPoints = mesh_.points();
393 
394  forAll(meshCells, celli)
395  {
396  bool found = false;
397 
398  for (const label facei : meshCells[celli])
399  {
400  for (const label fpi : meshFaces[facei])
401  {
402  if (bb.contains(meshPoints[fpi]))
403  {
404  found = true;
405  break;
406  }
407  }
408 
409  if (found)
410  {
411  treeCellIDs.append(celli);
412  break;
413  }
414  }
415  }
416 
418  (
419  treeDataCell(true, mesh_, std::move(treeCellIDs), polyMesh::CELL_TETS),
420  bb,
421  10,
422  100,
423  10
424  );
425 
426  cellIds_.setSize(points_.size(), -1);
427  pointMask_.setSize(points_.size(), false);
428 
429  // Kick the tet base points calculation to avoid parallel comms later
430  (void)mesh_.tetBasePtIs();
431 
432  const auto& cellLabels = tree.shapes().cellLabels();
433 
434  forAll(points_, pointi)
435  {
436  const vector& pos = points_[pointi];
437 
438 // label meshCelli = mesh_.findCell(pos);
439  label treeCelli = tree.findInside(pos);
440 
441  label proci = treeCelli >= 0 ? Pstream::myProcNo() : -1;
442 
443  reduce(proci, maxOp<label>());
444 
445  pointMask_[pointi] = treeCelli != -1;
446 
447  if (proci >= 0)
448  {
449  cellIds_[pointi] =
450  proci == Pstream::myProcNo()
451  ? cellLabels[treeCelli]
452  : -1;
453  }
454  else
455  {
456  if (errorOnPointNotFound_)
457  {
459  << "Position " << pos << " not found in mesh"
460  << abort(FatalError);
461  }
462  else
463  {
464  DebugInfo
465  << "Position " << pos << " not found in mesh"
466  << endl;
467  }
468  }
469  }
470 
472  Pstream::listCombineScatter(pointMask_);
473 }
474 
475 
477 (
478  const scalarField& field
479 ) const
480 {
481  if (field.size() != points_.size())
482  {
484  << "Inconsistent field sizes: input:" << field.size()
485  << " points:" << points_.size()
486  << abort(FatalError);
487  }
488 
489  scalar sumArea = 0;
490  scalar sumFieldArea = 0;
491  forAll(faces_, facei)
492  {
493  const face& f = faces_[facei];
494 
495  bool valid = true;
496  scalar faceValue = 0;
497  for (const label pti : f)
498  {
499  // Exclude contributions where sample cell for point was not found
500  if (!pointMask_[pti])
501  {
502  valid = false;
503  break;
504  }
505  faceValue += field[pti];
506  }
507 
508  if (valid)
509  {
510  scalar area = f.mag(points_);
511  sumArea += area;
512  sumFieldArea += faceValue/f.size()*area;
513  }
514  }
515 
516  return sumFieldArea/(sumArea + ROOTVSMALL);
517 }
518 
519 
521 (
522  const vectorField& U,
523  const vectorField& Ur,
524  const scalar URef
525 )
526 {
527  // Write surface
528  if (!surfaceWriterPtr_)
529  {
530  return;
531  }
532 
533 
534  // Time-aware, with time spliced into the output path
535  surfaceWriterPtr_->isPointData(true);
536  surfaceWriterPtr_->beginTime(time_);
537  surfaceWriterPtr_->open
538  (
539  points_,
540  faces_,
541  (baseFileDir() / name() / "surfaces" / "disk"),
542  false // serial - already merged
543  );
544  surfaceWriterPtr_->nFields(4); // Legacy VTK
545  surfaceWriterPtr_->write("U", U);
546  surfaceWriterPtr_->write("Ur", Ur);
547  surfaceWriterPtr_->write("UNorm", U/URef);
548  surfaceWriterPtr_->write("UrNorm", Ur/URef);
549  surfaceWriterPtr_->endTime();
550  surfaceWriterPtr_->clear();
551 }
552 
553 
555 {
556  if (!writePropellerPerformance_)
557  {
558  return;
559  }
560 
561  // Update n_
562  setRotationalSpeed();
563 
564  const vector sumForce(sum(force_[0]) + sum(force_[1]) + sum(force_[2]));
565  const vector sumMoment(sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]));
566 
567  const scalar diameter = 2*radius_;
568  const scalar URef = URefPtr_->value(time_.timeOutputValue());
569  const scalar j = -URef/mag(n_ + ROOTVSMALL)/diameter;
570  const scalar denom = rhoRef_*sqr(n_)*pow4(diameter);
571  const scalar kt = (sumForce & coordSysPtr_->e3())/denom;
572  const scalar kq =
573  -sign(n_)*(sumMoment & coordSysPtr_->e3())/(denom*diameter);
574  const scalar etaO = kt*j/(kq*constant::mathematical::twoPi + ROOTVSMALL);
575 
576  if (writeToFile())
577  {
578  auto& os = propellerPerformanceFilePtr_();
579 
580  writeCurrentTime(os);
581  os << tab << n_
582  << tab << URef
583  << tab << j
584  << tab << kt
585  << tab << 10*kq
586  << tab << etaO
587  << nl;
588 
589  os.flush();
590  }
591 
592  Log << type() << " " << name() << " output:" << nl
593  << " Revolutions per second, n : " << n_ << nl
594  << " Reference velocity, URef : " << URef << nl
595  << " Advance coefficient, J : " << j << nl
596  << " Thrust coefficient, Kt : " << kt << nl
597  << " Torque coefficient, 10*Kq : " << 10*kq << nl
598  << " Efficiency, etaO : " << etaO << nl
599  << nl;
600 
601 
602  // Write state/results information
603  setResult("n", n_);
604  setResult("URef", URef);
605  setResult("Kt", kt);
606  setResult("Kq", kq);
607  setResult("J", j);
608  setResult("etaO", etaO);
609 }
610 
611 
613 (
614  const vectorField& U,
615  const scalar URef
616 )
617 {
618  if (!Pstream::master()) return;
619 
620  // Velocity
621  auto& os = wakeFilePtr_();
622 
623  const pointField propPoints(coordSysPtr_->localPosition(points_));
624  const label offset =
625  mag(propPoints[1][0] - propPoints[0][0]) < SMALL ? 0 : 1;
626  const scalar rMax = propPoints.last()[0];
627 
628  const scalar UzMean = meanSampleDiskField(U.component(2));
629 
630  writeHeaderValue(os, "Time", time_.timeOutputValue());
631  writeHeaderValue(os, "Reference velocity", URef);
632  writeHeaderValue(os, "Direction", coordSysPtr_->e3());
633  writeHeaderValue(os, "Wake", 1 - UzMean/URef);
634  writeHeader(os, "");
635  writeCommented(os, "r/R");
636  writeTabbed(os, "alpha");
637  writeTabbed(os, "(x y z)");
638  writeTabbed(os, "(Ur Utheta Uz)");
639  os << nl;
640 
641  for (label thetai = 0; thetai < nTheta_; ++thetai)
642  {
643  const scalar deg = 360*thetai/scalar(nTheta_);
644 
645  for (label radiusi = 0; radiusi <= nRadial_; ++radiusi)
646  {
647  label pointi = radiusi*nTheta_ + thetai + offset;
648 
649  if (radiusi == 0 && offset == 1)
650  {
651  // Only a single point at the centre - repeat for all thetai
652  pointi = 0;
653  }
654 
655  if (pointMask_[pointi])
656  {
657  const scalar rR = propPoints[radiusi*nTheta_][0]/rMax;
658 
659  os << rR << tab << deg << tab
660  << points_[pointi] << tab << U[pointi] << nl;
661  }
662  }
663  }
664 
665  writeBreak(os);
666 
667  os << endl;
668 }
669 
670 
672 (
673  const vectorField& U,
674  const scalar URef
675 )
676 {
677  if (!Pstream::master()) return;
678 
679  // Alternative common format - axial wake component
680  auto& os = axialWakeFilePtr_();
681 
682  const pointField propPoints(coordSysPtr_->localPosition(points_));
683  const label offset =
684  mag(propPoints[1][0] - propPoints[0][0]) < SMALL ? 0 : 1;
685  const scalar rMax = propPoints.last()[0];
686 
687  writeHeaderValue(os, "Time", time_.timeOutputValue());
688 
689  os << "# angle";
690  for (label radiusi = 0; radiusi <= nRadial_; ++radiusi)
691  {
692  label pointi = radiusi*nTheta_;
693  scalar r = propPoints[pointi][0];
694  os << tab << "r/R=" << r/rMax;
695  }
696  os << nl;
697 
698  for (label thetai = 0; thetai < nTheta_; ++thetai)
699  {
700  os << 360*thetai/scalar(nTheta_);
701 
702  for (label radiusi = 0; radiusi <= nRadial_; ++radiusi)
703  {
704  label pointi = radiusi*nTheta_ + thetai + offset;
705 
706  if (radiusi == 0 && offset == 1)
707  {
708  // Only a single point at the centre - repeat for all thetai
709  pointi = 0;
710  }
711 
712  if (pointMask_[pointi])
713  {
714  os << tab << 1 - U[pointi][2]/URef;
715  }
716  else
717  {
718  os << tab << "undefined";
719  }
720  }
721 
722  os << nl;
723  }
724 
725  writeBreak(os);
726 
727  os << endl;
728 }
729 
730 
732 {
733  if (!writeWakeFields_)
734  {
735  return;
736  }
737 
738  // Normalised velocity
739  const vectorField UDisk(interpolate(U(), vector::uniform(nanValue_))());
740  const vectorField UrDisk(coordSysPtr_->localVector(UDisk));
741 
742  // Surface field data
743  writeSampleDiskSurface(UDisk, UrDisk, URef);
744 
745  // Write wake text files
746  writeWake(UrDisk, URef);
747  writeAxialWake(UrDisk, URef);
748 }
749 
750 
751 template<class Type>
753 (
755  const Type& defaultValue
756 ) const
757 {
758  auto tfield = tmp<Field<Type>>::New(points_.size(), defaultValue);
759  auto& field = tfield.ref();
760 
761  autoPtr<interpolation<Type>> interpolator
762  (
763  interpolation<Type>::New(interpolationScheme_, psi)
764  );
765 
766  forAll(points_, pointi)
767  {
768  const label celli = cellIds_[pointi];
769 
770  if (cellIds_[pointi] != -1)
771  {
772  const point& position = points_[pointi];
773  field[pointi] = interpolator().interpolate(position, celli);
774  }
775  }
776 
779 
780  return tfield;
781 }
782 
783 
784 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
785 
787 (
788  const word& name,
789  const Time& runTime,
790  const dictionary& dict,
791  bool readFields
792 )
793 :
794  forces(name, runTime, dict, false),
795  radius_(0),
796  URefPtr_(nullptr),
797  rotationMode_(rotationMode::SPECIFIED),
798  MRFName_(),
799  n_(0),
800  writePropellerPerformance_(true),
801  propellerPerformanceFilePtr_(nullptr),
802  writeWakeFields_(true),
803  surfaceWriterPtr_(nullptr),
804  nTheta_(0),
805  nRadial_(0),
806  points_(),
807  errorOnPointNotFound_(false),
808  faces_(),
809  cellIds_(),
810  pointMask_(),
811  interpolationScheme_("cell"),
812  wakeFilePtr_(nullptr),
813  axialWakeFilePtr_(nullptr),
814  nanValue_(pTraits<scalar>::min)
815 {
816  if (readFields)
817  {
818  read(dict);
819  Log << endl;
820  }
821 }
822 
823 
825 (
826  const word& name,
827  const objectRegistry& obr,
828  const dictionary& dict,
829  bool readFields
830 )
831 :
832  propellerInfo(name, obr.time(), dict, false)
833 {}
834 
835 
836 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
837 
839 {
840  if (forces::read(dict))
841  {
842  radius_ = dict.getScalar("radius");
843  URefPtr_.reset(Function1<scalar>::New("URef", dict, &mesh_));
844  rotationMode_ = rotationModeNames_.get("rotationMode", dict);
845 
846  // Must be set before setting the surface
847  setCoordinateSystem(dict);
848 
849  writePropellerPerformance_ =
850  dict.get<bool>("writePropellerPerformance");
851 
852  writeWakeFields_ = dict.get<bool>("writeWakeFields");
853  if (writeWakeFields_)
854  {
855  setSampleDiskSurface(dict);
856 
857  dict.readIfPresent("interpolationScheme", interpolationScheme_);
858 
859  dict.readIfPresent("nanValue", nanValue_);
860  }
861 
862  return true;
863  }
864 
865  return false;
866 }
867 
868 
870 {
871  calcForcesMoment();
872 
873  createFiles();
874 
875  if (writeWakeFields_)
876  {
877  // Only setting mean axial velocity result during execute
878  // - wake fields are 'heavy' and controlled separately using the
879  // writeControl
880  const vectorField
881  UDisk
882  (
883  coordSysPtr_->localVector
884  (
886  (
887  U(),
888  vector::uniform(nanValue_)
889  )()
890  )
891  );
892  const scalar UzMean = meanSampleDiskField(UDisk.component(2));
893 
894  setResult("UzMean", UzMean);
895  }
896 
897  writePropellerPerformance();
898 
899  return true;
900 }
901 
902 
904 {
905  const scalar URef = URefPtr_->value(time_.timeOutputValue());
906  writeWakeFields(URef);
907 
908  return true;
909 }
910 
911 
913 {
914  updateSampleDiskCells();
915 }
916 
917 
919 {
920  updateSampleDiskCells();
921 }
922 
923 
924 // ************************************************************************* //
Foam::coordinateSystem::globalPosition
point globalPosition(const point &local) const
From local coordinate position to global (cartesian) position.
Definition: coordinateSystem.H:595
Foam::functionObjects::propellerInfo::movePoints
void movePoints(const polyMesh &mesh)
Update for changes of mesh.
Definition: propellerInfo.C:918
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::maxOp
Definition: ops.H:223
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::uniform
static Vector< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:164
Foam::dictionary::getScalar
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get< scalar >(const word&, keyType::option)
Definition: dictionary.H:1512
Foam::functionObjects::forces
Calculates the forces and moments by integrating the pressure and skin-friction forces over a given l...
Definition: forces.H:236
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< label >
Foam::boundBox::inflate
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:175
Function1.H
Foam::functionObjects::propellerInfo::createFiles
void createFiles()
Create output files.
Definition: propellerInfo.C:175
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
propellerInfo.H
Foam::functionObjects::propellerInfo::writePropellerPerformance
void writePropellerPerformance()
Write the wake fields.
Definition: propellerInfo.C:554
Foam::functionObjects::propellerInfo::execute
virtual bool execute()
Execute, currently does nothing.
Definition: propellerInfo.C:869
Foam::functionObjects::propellerInfo::MRFName_
word MRFName_
Name of MRF zone (if applicable)
Definition: propellerInfo.H:332
Foam::functionObjects::propellerInfo::rotationMode::MRF
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::indexedOctree::shapes
const Type & shapes() const
Reference to shape.
Definition: indexedOctree.H:444
Foam::orEqOp
Definition: ops.H:86
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
interpolation.H
Foam::writeHeader
static void writeHeader(Ostream &os, const word &fieldName)
Definition: rawSurfaceWriterImpl.C:66
surfaceWriter.H
Foam::functionObjects::propellerInfo::rotationMode::SPECIFIED
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::DynamicList::clear
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
Foam::functionObjects::propellerInfo::writeAxialWake
void writeAxialWake(const vectorField &U, const scalar URef)
Write the axial wake text file.
Definition: propellerInfo.C:672
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::functionObjects::propellerInfo::propellerInfo
propellerInfo(const propellerInfo &)=delete
No copy construct.
Foam::functionObjects::propellerInfo::meanSampleDiskField
scalar meanSampleDiskField(const scalarField &field) const
Return the area average of a field.
Definition: propellerInfo.C:477
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::functionObjects::propellerInfo::writeWakeFields
void writeWakeFields(const scalar URef)
Write the wake fields.
Definition: propellerInfo.C:731
Foam::coordinateSystem::origin
virtual const point & origin() const
Return origin.
Definition: coordinateSystem.H:469
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::functionObjects::propellerInfo::write
virtual bool write()
Write the forces.
Definition: propellerInfo.C:903
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
Foam::treeDataCell
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
Definition: treeDataCell.H:56
Foam::Field< vector >
Foam::functionObjects::propellerInfo::U
const volVectorField & U() const
Return the velocity field.
Definition: propellerInfo.C:212
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::functionObjects::propellerInfo::updateSampleDiskCells
void updateSampleDiskCells()
Set the sample cells corresponding to the sample points.
Definition: propellerInfo.C:379
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
Foam::functionObjects::propellerInfo
Calculates propeller performance and wake field properties.
Definition: propellerInfo.H:302
field
rDeltaTY field()
Foam::functionObjects::forces::read
virtual bool read(const dictionary &)
Read the forces data.
Definition: forces.C:799
Foam::functionObjects::propellerInfo::rotationModeNames_
static const Enum< rotationMode > rotationModeNames_
Definition: propellerInfo.H:315
cylindricalCS.H
Foam::indexedOctree::findInside
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
Definition: indexedOctree.C:2629
Foam::coordSystem::cylindrical
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
Definition: cylindricalCS.H:71
Foam::constant::mathematical::twoPi
constexpr scalar twoPi(2 *M_PI)
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
treeDataCell.H
Foam::functionObjects::propellerInfo::n_
scalar n_
Propeller speed (revolutions per second)
Definition: propellerInfo.H:335
Foam::functionObjects::propellerInfo::writeSampleDiskSurface
void writeSampleDiskSurface(const vectorField &U, const vectorField &Ur, const scalar URef)
Write the sample surface.
Definition: propellerInfo.C:521
Foam::maxEqOp
Definition: ops.H:80
Foam::functionObjects::propellerInfo::setCoordinateSystem
void setCoordinateSystem(const dictionary &dict)
Set the coordinate system.
Definition: propellerInfo.C:61
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
IOMRFZoneList.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::polyMesh::CELL_TETS
Definition: polyMesh.H:109
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::treeBoundBox::contains
bool contains(const vector &dir, const point &) const
Contains point (inside or on edge) and moving in direction.
Definition: treeBoundBox.C:320
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::objectRegistry::cfindObject
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:390
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::propellerInfo::rotationMode_
rotationMode rotationMode_
Rotation mode.
Definition: propellerInfo.H:329
f
labelList f(nPoints)
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::surfaceWriter::New
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:64
Foam::Vector< scalar >
Foam::List< face >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::functionObjects::propellerInfo::writeWake
void writeWake(const vectorField &U, const scalar URef)
Write the wake text file.
Definition: propellerInfo.C:613
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::IOMRFZoneList
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries,...
Definition: IOMRFZoneList.H:71
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::functionObjects::propellerInfo::setSampleDiskSurface
void setSampleDiskSurface(const dictionary &dict)
Set the sample surface based on dictionary settings.
Definition: propellerInfo.C:333
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:155
Foam::functionObjects::propellerInfo::setSampleDiskGeometry
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.
Definition: propellerInfo.C:230
Foam::functionObjects::propellerInfo::UpdateMesh
void UpdateMesh(const mapPolyMesh &mpm)
Definition: propellerInfo.C:912
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::functionObjects::propellerInfo::read
virtual bool read(const dictionary &)
Read the forces data.
Definition: propellerInfo.C:838
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::functionObjects::propellerInfo::interpolate
tmp< Field< Type > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &psi, const Type &defaultValue) const
Interpolate from the mesh onto the sample surface.
Foam::GeometricField< vector, fvPatchField, volMesh >
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::functionObjects::propellerInfo::setRotationalSpeed
void setRotationalSpeed()
Set the rotational speed.
Definition: propellerInfo.C:137
Foam::coordinateSystem
Base class for coordinate system specification, the default coordinate system type is cartesian .
Definition: coordinateSystem.H:132
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177