lumpedPointMovement.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) 2016-2023 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 "lumpedPointMovement.H"
29 #include "lumpedPointIOMovement.H"
30 #include "Fstream.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "PtrMap.H"
34 #include "triFace.H"
35 #include "labelPair.H"
36 #include "indexedOctree.H"
37 #include "treeDataPoint.H"
38 #include "pointIndexHit.H"
39 #include "pointPatch.H"
40 #include "PstreamReduceOps.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
45 (
46  ::Foam::debug::debugSwitch("lumpedPointMovement", 0)
47 );
48 
49 
50 const Foam::word
51 Foam::lumpedPointMovement::canonicalName("lumpedPointMovement");
52 
53 
54 const Foam::Enum
55 <
57 >
59 ({
60  { outputFormatType::PLAIN, "plain" },
61  { outputFormatType::DICTIONARY, "dictionary" },
62 });
63 
64 
65 const Foam::Enum
66 <
68 >
70 ({
71  { scalingType::LENGTH, "length" },
72  { scalingType::FORCE, "force" },
73  { scalingType::MOMENT, "moment" },
74 });
75 
76 
77 
78 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
79 
80 namespace Foam
81 {
82 
84 //- Space-separated vector value (ASCII)
85 static inline Ostream& putPlain(Ostream& os, const vector& v)
86 {
87  return os << v.x() << ' ' << v.y() << ' ' << v.z();
88 }
89 
90 
92 //- Write list content with size, bracket, content, bracket one-per-line.
93 // This makes for consistent for parsing, regardless of the list length.
94 template <class T>
95 static void writeList(Ostream& os, const string& header, const UList<T>& list)
96 {
97  const label len = list.size();
98 
99  // Header string
100  os << header.c_str() << nl;
101 
102  // Write size and start delimiter
103  os << len << nl << token::BEGIN_LIST << nl;
104 
105  // Write contents
106  for (label i=0; i < len; ++i)
107  {
108  os << list[i] << nl;
109  }
110 
111  // Write end delimiter
112  os << token::END_LIST;
113  os.endEntry() << nl;
114 }
116 
117 } // End namespace Foam
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 :
124  origin_(Zero),
125  state0_(),
126  state_(),
127  originalIds_(),
128  controllers_(),
129  patchControls_(),
130  relax_(1),
131  ownerId_(-1),
132  forcesDict_(),
133  coupler_(),
134  inputName_("positions.in"),
135  outputName_("forces.out"),
136  logName_("movement.log"),
137  inputFormat_(lumpedPointState::inputFormatType::DICTIONARY),
138  outputFormat_(outputFormatType::DICTIONARY),
139  scaleInput_(-1),
140  scaleOutput_(-1),
141  calcFrequency_(1),
142  lastTrigger_(-1)
143 {}
144 
147 (
148  const dictionary& dict,
149  label ownerId
150 )
151 :
152  origin_(Zero),
153  state0_(),
154  state_(),
155  originalIds_(),
156  controllers_(),
157  patchControls_(),
158  relax_(1),
159  ownerId_(ownerId),
160  forcesDict_(),
161  coupler_(),
162  inputName_("positions.in"),
163  outputName_("forces.out"),
164  logName_("movement.log"),
165  scaleInput_(-1),
166  scaleOutput_(-1),
167  calcFrequency_(1),
168  lastTrigger_(-1)
169 {
170  readDict(dict);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
177 {
178  return (timeIndex >= lastTrigger_ + calcFrequency_);
179 }
181 
183 {
184  lastTrigger_ = timeIndex;
185 }
187 
189 {
190  origin_ = dict.getOrDefault<point>("origin", Zero);
191 
192  // Initial point locations (zero rotation)
193  auto tpoints0 = tmp<pointField>::New();
194  auto& points0 = tpoints0.ref();
195 
196  dict.readEntry("points", points0);
197  points0 += origin_;
198 
199  originalIds_.clear();
200  controllers_.clear();
201  patchControls_.clear();
202 
203 
204  // The FEA ids for the points (optional)
205  Map<label> pointIdMap;
206 
207  if (dict.readIfPresent("pointLabels", originalIds_) && originalIds_.size())
208  {
209  if (originalIds_.size() != points0.size())
210  {
212  << "Incorrect number of pointLabels. Had "
213  << originalIds_.size() << " for " << points0.size() << " points"
214  << nl << endl
215  << exit(FatalIOError);
216  }
217 
218  labelHashSet duplicates;
219 
220  label pointi = 0;
221 
222  for (const label id : originalIds_)
223  {
224  if (!pointIdMap.insert(id, pointi))
225  {
226  duplicates.set(id);
227  }
228 
229  ++pointi;
230  }
231 
232  if (!duplicates.empty())
233  {
235  << "Found duplicate point ids "
236  << flatOutput(duplicates.sortedToc()) << nl << endl
237  << exit(FatalIOError);
238  }
239  }
240 
241  const dictionary* dictptr = dict.findDict("controllers");
242  if (dictptr)
243  {
244  controllers_ = HashPtrTable<lumpedPointController>(*dictptr);
245 
246  // Verify the input
247  forAllIters(controllers_, iter)
248  {
249  (*iter)->remapPointLabels(points0.size(), pointIdMap);
250  }
251  }
252  else
253  {
254  // Add single global controller
255  // Warning?
256 
257  controllers_.clear();
258  }
259 
260  relax_ = dict.getOrDefault<scalar>("relax", 1);
261  relax_ = clamp(relax_, zero_one{});
262 
263  forcesDict_.merge(dict.subOrEmptyDict("forces"));
264 
265  const dictionary& commDict = dict.subDict("communication");
266  coupler_.readDict(commDict);
267 
268  calcFrequency_ = commDict.getOrDefault<label>("calcFrequency", 1);
269 
270  // Leave trigger intact
271 
272  commDict.readEntry("inputName", inputName_);
273  commDict.readEntry("outputName", outputName_);
274  commDict.readIfPresent("logName", logName_);
275 
276  inputFormat_ = lumpedPointState::formatNames.get("inputFormat", commDict);
277  outputFormat_ = formatNames.get("outputFormat", commDict);
278 
279  scaleInput_ = -1;
280  scaleOutput_ = -1;
281 
282  const dictionary* scaleDict = nullptr;
283 
284  if ((scaleDict = commDict.findDict("scaleInput")) != nullptr)
285  {
286  for (int i=0; i < scaleInput_.size(); ++i)
287  {
288  const word& key = scalingNames[scalingType(i)];
289 
290  if
291  (
292  scaleDict->readIfPresent(key, scaleInput_[i])
293  && scaleInput_[i] > 0
294  )
295  {
296  Info<<"Using input " << key << " multiplier: "
297  << scaleInput_[i] << nl;
298  }
299  }
300  }
301 
302  if ((scaleDict = commDict.findDict("scaleOutput")) != nullptr)
303  {
304  for (int i=0; i < scaleOutput_.size(); ++i)
305  {
306  const word& key = scalingNames[scalingType(i)];
307 
308  if
309  (
310  scaleDict->readIfPresent(key, scaleOutput_[i])
311  && scaleOutput_[i] > 0
312  )
313  {
314  Info<<"Using output " << key << " multiplier: "
315  << scaleOutput_[i] << nl;
316  }
317  }
318  }
319 
320  state0_ = lumpedPointState
321  (
322  tpoints0,
323  quaternion::eulerOrderNames.getOrDefault
324  (
325  "rotationOrder",
326  dict,
327  quaternion::eulerOrder::ZXZ
328  ),
329  dict.getOrDefault("degrees", false)
330  );
331 
332  state0_.scalePoints(scaleInput_[scalingType::LENGTH]);
333 
334  state_ = state0_;
335 }
336 
339 (
340  const polyPatch& pp
341 ) const
342 {
343  return hasPatchControl(pp.index());
344 }
345 
348 (
349  const pointPatch& fpatch
350 ) const
351 {
352  return hasInterpolator(fpatch.index());
353 }
354 
357 (
358  const polyPatch& pp
359 ) const
360 {
361  const auto ctrlIter = patchControls_.cfind(pp.index());
362 
363  if (!ctrlIter.good())
364  {
366  << "No controllers for patch " << pp.name()
367  << exit(FatalError);
368  }
369 
370  const patchControl& ctrl = *ctrlIter;
371 
372  for (const word& ctrlName : ctrl.names_)
373  {
374  const auto iter = controllers_.cfind(ctrlName);
375 
376  if (!iter.good())
377  {
379  << "No controller: " << ctrlName << nl
380  << " For patch " << pp.name()
381  << exit(FatalError);
382  }
383  }
384 }
385 
388 (
389  const polyPatch& pp,
390  const wordList& ctrlNames,
391  const pointField& points0
392 )
393 {
394  // Reference mass centres
395  const pointField& lumpedCentres0 = state0().points();
396 
397  const label patchIndex = pp.index();
398 
399  // Info<<"Add patch control for patch " << patchIndex << " "
400  // << pp.name() << nl;
401 
402  patchControl& ctrl = patchControls_(patchIndex);
403  ctrl.names_ = ctrlNames;
404 
405  labelList& faceToPoint = ctrl.faceToPoint_;
406  faceToPoint.resize(pp.size(), -1);
407 
408  checkPatchControl(pp);
409 
410  const faceList& faces = pp.boundaryMesh().mesh().faces();
411 
412  // Subset of points to search (if specified)
413  labelHashSet subsetPointIds;
414 
415  for (const word& ctrlName : ctrl.names_)
416  {
417  const auto iter = controllers_.cfind(ctrlName);
418 
419  if (!iter.good())
420  {
422  << "No controller: " << ctrlName << nl
423  << exit(FatalError);
424  }
425 
426  const labelList& pointLabels = (*iter)->pointLabels();
427 
428  subsetPointIds.insert(pointLabels);
429  }
430 
431  if (ctrl.names_.size() && subsetPointIds.empty())
432  {
434  << "Controllers specified, but without any points" << nl
435  << exit(FatalError);
436  }
437 
438 
439  treeBoundBox bb(lumpedCentres0);
440  bb.inflate(0.01);
441 
443  (
445  (
446  lumpedCentres0,
447  subsetPointIds.sortedToc(),
448  subsetPointIds.size() // subset is optional/required
449  ),
450  bb, // overall search domain
451  8, // maxLevel
452  10, // leafsize
453  3.0 // duplicity
454  );
455  const auto& treePoints = ppTree.shapes();
456 
457  const scalar searchDistSqr(sqr(GREAT));
458 
459  const label patchStart = pp.start();
460  forAll(pp, patchFacei)
461  {
462  const point fc(faces[patchStart + patchFacei].centre(points0));
463 
464  // Store the original pointId, not subset id
465  faceToPoint[patchFacei] =
466  treePoints.objectIndex
467  (
468  ppTree.findNearest(fc, searchDistSqr).index()
469  );
470  }
471 
472  if (debug)
473  {
474  Pout<<"Added face mapping for patch: " << patchIndex << endl;
475  }
476 }
477 
480 (
481  const pointPatch& fpatch,
482  const pointField& points0
483 )
484 {
485  // Reference mass centres
486  const pointField& lumpedCentres0 = state0().points();
487 
488  const label patchIndex = fpatch.index();
489 
490  patchControl& ctrl = patchControls_(patchIndex);
491 
492  List<lumpedPointInterpolator>& interpList = ctrl.interp_;
493  interpList.clear();
494 
495  // The connectivity, adjacency list
496  Map<labelHashSet> adjacency;
497 
498  // Subset of points to search (if specified)
499  labelHashSet subsetPointIds;
500 
501  for (const word& ctrlName : ctrl.names_)
502  {
503  const auto iter = controllers_.cfind(ctrlName);
504 
505  if (!iter.good())
506  {
508  << "No controller: " << ctrlName << nl
509  << exit(FatalError);
510  }
511 
512  const labelList& pointLabels = (*iter)->pointLabels();
513 
514  subsetPointIds.insert(pointLabels);
515 
516  // Adjacency lists
517  forAll(pointLabels, i)
518  {
519  const label curr = pointLabels[i];
520 
521  if (i)
522  {
523  const label prev = pointLabels[i-1];
524  adjacency(prev).insert(curr);
525  adjacency(curr).insert(prev);
526  }
527  else if (!adjacency.found(curr))
528  {
529  adjacency(curr).clear();
530  }
531  }
532  }
533 
534  if (ctrl.names_.empty())
535  {
536  // Adjacency lists
537  const label len = state0().size();
538 
539  for (label i=0; i < len; ++i)
540  {
541  const label curr = i;
542 
543  if (i)
544  {
545  const label prev = i-1;
546  adjacency(prev).insert(curr);
547  adjacency(curr).insert(prev);
548  }
549  else if (!adjacency.found(curr))
550  {
551  adjacency(curr).clear();
552  }
553  }
554  }
555 
556  if (ctrl.names_.size() && subsetPointIds.empty())
557  {
559  << "Controllers specified, but without any points" << nl
560  << exit(FatalError);
561  }
562 
563 
564  // Pairs defining connecting points as triangles
565  Map<labelPairList> adjacencyPairs;
566 
567  barycentric2D bary;
568 
569  {
570  // Pairs for the ends
572 
573  // Mag sin(angle) around the triangle point 0,
574  // used to sort the generated triangles according to the acute angle
575  DynamicList<scalar> acuteAngles;
576 
577  forAllConstIters(adjacency, iter)
578  {
579  const label nearest = iter.key();
580 
581  labelList neighbours(iter.val().sortedToc());
582 
583  const label len = neighbours.size();
584 
585  pairs.clear();
586  acuteAngles.clear();
587 
588  const point& nearPt = lumpedCentres0[nearest];
589 
590  for (label j=1; j < len; ++j)
591  {
592  for (label i=0; i < j; ++i)
593  {
594  labelPair neiPair(neighbours[i], neighbours[j]);
595 
596  triPointRef tri
597  (
598  nearPt,
599  lumpedCentres0[neiPair.first()],
600  lumpedCentres0[neiPair.second()]
601  );
602 
603  // Require non-degenerate triangles
604  if (tri.pointToBarycentric(tri.a(), bary) > SMALL)
605  {
606  // Triangle OK
607  pairs.append(neiPair);
608 
609  vector ab(normalised(tri.b() - tri.a()));
610  vector ac(normalised(tri.c() - tri.a()));
611 
612  // Angle between neighbouring edges
613  // Use negative cosine to map [0-180] -> [-1 .. +1]
614 
615  acuteAngles.append(-(ab & ac));
616  }
617  }
618  }
619 
620  if (pairs.size() > 1)
621  {
622  // Sort by acute angle
623  labelList order(sortedOrder(acuteAngles));
624  inplaceReorder(order, pairs);
625  }
626 
627  adjacencyPairs.insert(nearest, pairs);
628  }
629  }
630 
631  if (debug & 2)
632  {
633  Info<< "Adjacency table for patch: " << fpatch.name() << nl;
634 
635  for (const label own : adjacency.sortedToc())
636  {
637  Info<< own << " =>";
638  for (const label nei : adjacency[own].sortedToc())
639  {
640  Info<< ' ' << nei;
641  }
642 
643  if (originalIds_.size())
644  {
645  Info<< " # " << originalIds_[own] << " =>";
646  for (const label nei : adjacency[own].sortedToc())
647  {
648  Info<< ' ' << originalIds_[nei];
649  }
650  }
651 
652  Info<< " # tri " << flatOutput(adjacencyPairs[own]);
653  Info<< nl;
654  }
655  }
656 
657  treeBoundBox bb(lumpedCentres0);
658  bb.inflate(0.01);
659 
661  (
663  (
664  lumpedCentres0,
665  subsetPointIds.sortedToc(),
666  subsetPointIds.size() // subset is optional/required
667  ),
668  bb, // overall search domain
669  8, // maxLevel
670  10, // leafsize
671  3.0 // duplicity
672  );
673  const auto& treePoints = ppTree.shapes();
674 
675 
676  // Searching
677 
678  const scalar searchDistSqr(sqr(GREAT));
679 
680  const labelList& meshPoints = fpatch.meshPoints();
681 
682  interpList.resize(meshPoints.size());
683 
684  DynamicList<scalar> unsortedNeiWeightDist;
685  DynamicList<label> unsortedNeighbours;
686 
687  forAll(meshPoints, pointi)
688  {
689  const point& ptOnMesh = points0[meshPoints[pointi]];
690 
691  // Nearest (original) point id
692  const label nearest =
693  treePoints.objectIndex
694  (
695  ppTree.findNearest(ptOnMesh, searchDistSqr).index()
696  );
697 
698  interpList[pointi].nearest(nearest);
699 
700  if (nearest == -1)
701  {
702  // Should not really happen
703  continue;
704  }
705 
706  // Have the nearest lumped point, now find the next nearest
707  // but check that the direction is also correct.
708 
709  // OK: within the 0-1 bounds
710  // 1+----+0
711  // .
712  // .
713  // +pt
714 
715  const point& nearPt = lumpedCentres0[nearest];
716 
717  const linePointRef toMeshPt(nearPt, ptOnMesh);
718 
719  const labelPairList& adjacentPairs = adjacencyPairs[nearest];
720 
721  unsortedNeighbours = adjacency[nearest].toc();
722  unsortedNeiWeightDist.resize(unsortedNeighbours.size());
723 
724  forAll(unsortedNeighbours, nbri)
725  {
726  unsortedNeiWeightDist[nbri] =
727  magSqr(ptOnMesh - lumpedCentres0[unsortedNeighbours[nbri]]);
728  }
729 
730  // Sort by distance
731  labelList distOrder(sortedOrder(unsortedNeiWeightDist));
732 
733  label ngood = 0;
734 
735  // Recalculate distance as weighting
736  for (const label nbri : distOrder)
737  {
738  const label nextPointi = unsortedNeighbours[nbri];
739 
740  const point& nextPt = lumpedCentres0[nextPointi];
741 
742  const linePointRef toNextPt(nearPt, nextPt);
743 
744  const scalar weight =
745  (toMeshPt.vec() & toNextPt.unitVec()) / toNextPt.mag();
746 
747  if (weight < 0)
748  {
749  // Reject: wrong direction or other bad value
750  continue;
751  }
752  else
753  {
754  // Store weight
755  unsortedNeiWeightDist[nbri] = weight;
756 
757  // Retain good weight
758  distOrder[ngood] = nbri;
759  ++ngood;
760  }
761  }
762 
763  distOrder.resize(ngood);
764 
765  if (distOrder.size() < 1)
766  {
767  continue;
768  }
769 
770  UIndirectList<label> neighbours(unsortedNeighbours, distOrder);
771  UIndirectList<scalar> neiWeight(unsortedNeiWeightDist, distOrder);
772 
773  bool useFirst = true;
774 
775  if (neighbours.size() > 1 && adjacentPairs.size())
776  {
777  // Check for best two neighbours
778 
779  bitSet neiPointid;
780  neiPointid.set(neighbours);
781 
782  for (const labelPair& ends : adjacentPairs)
783  {
784  label nei1 = ends.first();
785  label nei2 = ends.second();
786 
787  if (!neiPointid.test(nei1) || !neiPointid.test(nei2))
788  {
789  // Reject, invalid combination for this point location
790  continue;
791  }
792  else if (neighbours.find(nei2) < neighbours.find(nei1))
793  {
794  // Order by distance, which is not really needed,
795  // but helps with diagnostics
796  std::swap(nei1, nei2);
797  }
798 
799 
800  triFace triF(nearest, nei1, nei2);
801 
802  if
803  (
804  triF.tri(lumpedCentres0).pointToBarycentric(ptOnMesh, bary)
805  > SMALL
806  && !bary.outside()
807  )
808  {
809  // Use barycentric weights
810  interpList[pointi].set(triF, bary);
811 
812  useFirst = false;
813  break;
814  }
815  }
816  }
817 
818  if (useFirst)
819  {
820  // Use geometrically closest neighbour
821  interpList[pointi].next(neighbours.first(), neiWeight.first());
822  }
823  }
824 }
825 
828 (
829  const polyMesh& pmesh
830 ) const
831 {
832  const label nLumpedPoints = state0().size();
833 
834  List<scalar> zoneAreas(nLumpedPoints, Zero);
835 
836  if (patchControls_.empty())
837  {
839  << "Attempted to calculate areas without setMapping()"
840  << nl << endl;
841  return zoneAreas;
842  }
843 
844  const polyBoundaryMesh& patches = pmesh.boundaryMesh();
845 
846  // fvMesh and has pressure field
847  if (isA<fvMesh>(pmesh))
848  {
849  const auto& mesh = refCast<const fvMesh>(pmesh);
850 
851  // Face areas (on patches)
852  const surfaceVectorField::Boundary& patchSf =
853  mesh.Sf().boundaryField();
854 
855  forAllConstIters(patchControls_, iter)
856  {
857  const label patchIndex = iter.key();
858  const patchControl& ctrl = iter.val();
859 
860  const labelList& faceToPoint = ctrl.faceToPoint_;
861 
862  const polyPatch& pp = patches[patchIndex];
863 
864  forAll(pp, patchFacei)
865  {
866  // Force from the patch-face into sum
867  const label pointIndex = faceToPoint[patchFacei];
868 
869  if (pointIndex < 0)
870  {
871  // Unmapped, for whatever reason?
872  continue;
873  }
874 
875  // Force from the patch-face into sum
876  zoneAreas[pointIndex] += mag(patchSf[patchIndex][patchFacei]);
877  }
878  }
879  }
880 
882 
883  return zoneAreas;
884 }
885 
888 (
889  const polyMesh& pmesh,
890  List<vector>& forces,
891  List<vector>& moments
892 ) const
893 {
894  const label nLumpedPoints = state0().size();
895 
896  forces.resize(nLumpedPoints);
897  moments.resize(nLumpedPoints);
898 
899  if (patchControls_.empty())
900  {
902  << "Attempted to calculate forces without setMapping()"
903  << nl << endl;
904 
905  forces.resize(nLumpedPoints, Zero);
906  moments.resize(nLumpedPoints, Zero);
907  return false;
908  }
909 
910  // Initialize with zero
911  forces = Zero;
912  moments = Zero;
913 
914  // Current mass centres
915  const pointField& lumpedCentres = state().points();
916 
917  const polyBoundaryMesh& patches = pmesh.boundaryMesh();
918 
919  const word pName(forcesDict_.getOrDefault<word>("p", "p"));
920  scalar pRef = forcesDict_.getOrDefault<scalar>("pRef", 0);
921  scalar rhoRef = forcesDict_.getOrDefault<scalar>("rhoRef", 1);
922 
923 
924  // Calculated force per patch - cache
925  PtrMap<vectorField> forceOnPatches;
926 
927  const volScalarField* pPtr = pmesh.findObject<volScalarField>(pName);
928 
929  // fvMesh and has pressure field
930  if (isA<fvMesh>(pmesh) && pPtr)
931  {
932  const auto& mesh = refCast<const fvMesh>(pmesh);
933  const volScalarField& p = *pPtr;
934 
935  // Face centres (on patches)
936  const surfaceVectorField::Boundary& patchCf = mesh.Cf().boundaryField();
937 
938  // Face areas (on patches)
939  const surfaceVectorField::Boundary& patchSf = mesh.Sf().boundaryField();
940 
941  // Pressure (on patches)
942  const volScalarField::Boundary& patchPress = p.boundaryField();
943 
944  // rhoRef if the pressure field is dynamic, i.e. p/rho otherwise 1
945  rhoRef = (p.dimensions() == dimPressure ? 1.0 : rhoRef);
946 
947  // Scale pRef by density for incompressible simulations
948  pRef /= rhoRef;
949 
950  forAllConstIters(patchControls_, iter)
951  {
952  const label patchIndex = iter.key();
953  const patchControl& ctrl = iter.val();
954 
955  const labelList& faceToPoint = ctrl.faceToPoint_;
956 
957  if (!forceOnPatches.found(patchIndex))
958  {
959  // Patch faces are +ve outwards,
960  // so the forces (exerted by fluid on solid)
961  // already have the correct sign
962  forceOnPatches.set
963  (
964  patchIndex,
965  (
966  rhoRef
967  * patchSf[patchIndex] * (patchPress[patchIndex] - pRef)
968  ).ptr()
969  );
970  }
971 
972  const vectorField& forceOnPatch = *forceOnPatches[patchIndex];
973 
974  const polyPatch& pp = patches[patchIndex];
975 
976  forAll(pp, patchFacei)
977  {
978  // Force from the patch-face into sum
979  const label pointIndex = faceToPoint[patchFacei];
980 
981  if (pointIndex < 0)
982  {
983  // Unmapped, for whatever reason?
984  continue;
985  }
986 
987  // Force from the patch-face into sum
988  forces[pointIndex] += forceOnPatch[patchFacei];
989 
990  // Effective torque arm:
991  // - translated into the lumped-points coordinate system
992  // prior to determining the distance
993  const vector lever
994  (
995  patchCf[patchIndex][patchFacei]
996  - lumpedCentres[pointIndex]
997  );
998 
999  // Moment from the patch-face into sum
1000  moments[pointIndex] += lever ^ forceOnPatch[patchFacei];
1001  }
1002  }
1003  }
1004  else
1005  {
1006  Info<<"No pressure field" << endl;
1007  }
1008 
1011 
1012  return true;
1013 }
1014 
1015 
1018 (
1019  const pointPatch& fpatch,
1020  const pointField& points0
1021 ) const
1022 {
1023  return pointsDisplacement(state(), fpatch, points0);
1024 }
1025 
1026 
1029 (
1030  const lumpedPointState& state,
1031  const pointPatch& fpatch,
1032  const pointField& points0
1033 ) const
1034 {
1035  const label patchIndex = fpatch.index();
1036 
1037  // Reference mass centres
1038  const pointField& lumpedCentres0 = state0().points();
1039 
1040  // Current mass centres
1041  const pointField& lumpedCentres = state.points();
1042 
1043  // The local-to-global transformation tensor
1044  const tensorField& localToGlobal = state.rotations();
1045 
1046  const labelList& meshPoints = fpatch.meshPoints();
1047 
1048  // Could also verify the sizes (state vs original)
1049 
1050  auto tdisp = tmp<pointField>::New(fpatch.size());
1051  auto& disp = tdisp.ref();
1052 
1053  // The interpolator
1054  const List<lumpedPointInterpolator>& interpList
1055  = patchControls_[patchIndex].interp_;
1056 
1057  forAll(meshPoints, pointi)
1058  {
1059  const lumpedPointInterpolator& interp = interpList[pointi];
1060 
1061  const point& p0 = points0[meshPoints[pointi]];
1062 
1063  const vector origin0 = interp.interpolate(lumpedCentres0);
1064  const vector origin = interp.interpolate(lumpedCentres);
1065  const tensor rotTensor = interp.interpolate(localToGlobal);
1066 
1067  disp[pointi] = (rotTensor & (p0 - origin0)) + origin - p0;
1068  }
1069 
1070  return tdisp;
1071 }
1072 
1073 
1076 (
1077  const lumpedPointState& state,
1078  const pointPatch& fpatch,
1079  const pointField& points0
1080 ) const
1081 {
1082  const label patchIndex = fpatch.index();
1083 
1084  // Reference mass centres
1085  const pointField& lumpedCentres0 = state0().points();
1086 
1087  // Current mass centres
1088  const pointField& lumpedCentres = state.points();
1089 
1090  // The local-to-global transformation tensor
1091  const tensorField& localToGlobal = state.rotations();
1092 
1093  const labelList& meshPoints = fpatch.meshPoints();
1094 
1095  // Could also verify the sizes (state vs original)
1096 
1097  auto tdisp = tmp<pointField>::New(fpatch.size());
1098  auto& disp = tdisp.ref();
1099 
1100  // The interpolator
1101  const List<lumpedPointInterpolator>& interpList =
1102  patchControls_[patchIndex].interp_;
1103 
1104  forAll(meshPoints, pointi)
1105  {
1106  const lumpedPointInterpolator& interp = interpList[pointi];
1107 
1108  const point& p0 = points0[meshPoints[pointi]];
1109 
1110  const vector origin0 = interp.interpolate(lumpedCentres0);
1111  const vector origin = interp.interpolate(lumpedCentres);
1112  const tensor rotTensor = interp.interpolate(localToGlobal);
1113 
1114  disp[pointi] = (rotTensor & (p0 - origin0)) + origin;
1115  }
1116 
1117  return tdisp;
1118 }
1120 
1122 {
1123  // os.writeEntry("axis", axis_);
1124  // os.writeEntry("locations", locations_);
1125 }
1127 
1129 {
1130  lumpedPointState prev = state_;
1131 
1132  const bool status = state_.readData
1133  (
1134  inputFormat_,
1135  coupler().resolveFile(inputName_),
1136  state0().rotationOrder(),
1137  state0().degrees()
1138  );
1139 
1140  scalePoints(state_);
1141 
1142  state_.relax(relax_, prev);
1143 
1144  return status;
1145 }
1146 
1149 (
1150  Ostream& os,
1151  const UList<vector>& forces,
1152  const UList<vector>& moments,
1153  const outputFormatType fmt,
1154  const Tuple2<scalar, scalar>* timesWritten
1155 ) const
1156 {
1157  const bool writeMoments = (moments.size() == forces.size());
1158 
1159  if (fmt == outputFormatType::PLAIN)
1160  {
1161  os <<"########" << nl;
1162  if (timesWritten)
1163  {
1164  os << "# Time value=" << timesWritten->first() << nl
1165  << "# Time prev=" << timesWritten->second() << nl;
1166  }
1167  os << "# size=" << this->size() << nl
1168  << "# columns (points) (forces)";
1169 
1170  if (writeMoments)
1171  {
1172  os << " (moments)";
1173  }
1174 
1175  os << nl;
1176 
1177  bool report = false;
1178  scalar scaleLength = scaleOutput_[scalingType::LENGTH];
1179  scalar scaleForce = scaleOutput_[scalingType::FORCE];
1180  scalar scaleMoment = scaleOutput_[scalingType::MOMENT];
1181 
1182  if (scaleLength > 0)
1183  {
1184  report = true;
1185  }
1186  else
1187  {
1188  scaleLength = 1.0;
1189  }
1190 
1191  if (scaleForce > 0)
1192  {
1193  report = true;
1194  }
1195  else
1196  {
1197  scaleForce = 1.0;
1198  }
1199 
1200  if (writeMoments)
1201  {
1202  if (scaleMoment > 0)
1203  {
1204  report = true;
1205  }
1206  else
1207  {
1208  scaleMoment = 1.0;
1209  }
1210  }
1211 
1212  if (report)
1213  {
1214  os <<"# scaling points=" << scaleLength
1215  <<" forces=" << scaleForce;
1216 
1217  if (writeMoments)
1218  {
1219  os <<" moments=" << scaleMoment;
1220  }
1221 
1222  os << nl;
1223  }
1224 
1225  os <<"########" << nl;
1226 
1227  forAll(state0().points(), i)
1228  {
1229  const point& pt = state0().points()[i];
1230 
1231  putPlain(os, scaleLength * pt) << ' ';
1232 
1233  if (i < forces.size())
1234  {
1235  const vector val(scaleForce * forces[i]);
1236  putPlain(os, val);
1237  }
1238  else
1239  {
1240  putPlain(os, vector::zero);
1241  }
1242 
1243  if (writeMoments)
1244  {
1245  os << ' ';
1246  if (i < moments.size())
1247  {
1248  const vector val(scaleMoment * moments[i]);
1249  putPlain(os, val);
1250  }
1251  else
1252  {
1253  putPlain(os, vector::zero);
1254  }
1255  }
1256 
1257  os << nl;
1258  }
1259  }
1260  else
1261  {
1262  // Make it easier for external programs to parse
1263  // - exclude the usual OpenFOAM 'FoamFile' header
1264  // - ensure lists have consistent format
1265 
1266  os <<"////////" << nl;
1267  if (timesWritten)
1268  {
1269  os.writeEntry("time", timesWritten->first());
1270  os.writeEntry("prevTime", timesWritten->second());
1271  }
1272  os << nl;
1273 
1274  writeList(os, "points", state0().points());
1275  writeList(os, "forces", forces);
1276 
1277  if (writeMoments)
1278  {
1279  writeList(os, "moments", moments);
1280  }
1281  }
1282 
1283  return true;
1284 }
1285 
1288 (
1289  const UList<vector>& forces,
1290  const UList<vector>& moments,
1291  const Tuple2<scalar, scalar>* timesWritten
1292 ) const
1293 {
1294  if (!Pstream::master())
1295  {
1296  return false;
1297  }
1298 
1299  // Regular output
1300  {
1301  OFstream os
1302  (
1303  coupler().resolveFile(outputName_)
1304  );
1305 
1306  writeData(os, forces, moments, outputFormat_, timesWritten);
1307  }
1308 
1309  // Log output
1310  {
1311  OFstream os
1312  (
1313  coupler().resolveFile(logName_),
1314  IOstreamOption(),
1316  );
1317 
1318  writeData(os, forces, moments, outputFormatType::PLAIN, timesWritten);
1319  }
1320 
1321  return true;
1322 }
1323 
1324 
1325 // ************************************************************************* //
Foam::surfaceFields.
EnumType get(const word &enumName) const
The enumeration corresponding to the given name.
Definition: Enum.C:68
Point unitVec() const
Return the unit vector (start-to-end)
Definition: lineI.H:122
const T & first() const noexcept
Access the first element.
Definition: Pair.H:137
bool set(const Key &key)
Same as insert (no value to overwrite)
Definition: HashSet.H:240
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const surfaceVectorField & Sf() const
Return cell face area vectors.
A triangle primitive used to calculate face normals and swept volumes. Uses referred points...
Definition: triangle.H:217
A line primitive.
Definition: line.H:52
List< scalar > areas(const polyMesh &pmesh) const
The areas for each pressure-zone.
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:504
Inter-processor communication reduction functions.
Holds (reference to) pointField. Encapsulation of data needed for octree searches. Used for searching for nearest point. No bounding boxes around points. Only overlaps and calcNearest are implemented, rest makes little sense.
Definition: treeDataPoint.H:58
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
tmp< pointField > pointsPosition(const lumpedPointState &state, const pointPatch &fpatch, const pointField &points0) const
The points absolute position according to specified state.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
barycentric2D pointToBarycentric(const point &pt) const
Calculate the barycentric coordinates from the given point.
Definition: triangleI.H:424
labelList pointLabels(nPoints, -1)
static void writeData(Ostream &os, const Type &val)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
scalar mag() const
The magnitude (length) of the line.
Definition: lineI.H:96
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
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
void setInterpolator(const pointPatch &fpatch, const pointField &points0)
Check if patch control exists for specified patch.
triPointRef tri(const UList< point > &points) const
Return the triangle.
Definition: triFaceI.H:167
bool found(const label &key) const
Same as contains()
Definition: HashTable.H:1354
label find(const T &val, label pos=0, label len=-1) const
Find index of the first occurrence of the value.
const Type & shapes() const noexcept
Reference to shape.
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
bool outside() const
True if any coordinate is negative.
void inflate(const scalar factor)
Expand box by factor*mag(span) in all dimensions.
Definition: boundBoxI.H:381
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Output to file stream, using an OSstream.
Definition: OFstream.H:49
static const Enum< outputFormatType > formatNames
Names for the output format types.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
int debugSwitch(const char *name, const int deflt=0)
Lookup debug switch or add default value.
Definition: debug.C:222
void resize(const label len)
Alter addressable list size, allocating new space if required while recovering old content...
Definition: DynamicListI.H:353
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A HashTable of pointers to objects of type <T> with a label key.
Definition: HashTableFwd.H:42
T interpolate(const UList< T > &input) const
Linear interpolated value between nearest and next locations.
Begin list [isseparator].
Definition: token.H:161
void couplingCompleted(const label timeIndex) const
Register that coupling is completed at this calcFrequency.
void readDict(const dictionary &dict)
Update settings from dictionary.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
T & first()
Access first element of the list, position [0].
A simple container for options an IOstream can normally have.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers...
Definition: HashPtrTable.H:51
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
bool writeData(Ostream &os, const UList< vector > &forces, const UList< vector > &moments, const outputFormatType fmt=outputFormatType::PLAIN, const Tuple2< scalar, scalar > *timesWritten=nullptr) const
Write points, forces, moments. Only call from the master process.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
void setPatchControl(const polyPatch &pp, const wordList &ctrlNames, const pointField &points0)
Define pressure-zones mapping for faces in the specified patches.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
static const Enum< eulerOrder > eulerOrderNames
The names for Euler-angle and Tait-Bryan angles, including "rollPitchYaw" and "yawPitchRoll" aliases...
Definition: quaternion.H:132
static const Enum< scalingType > scalingNames
Names for the scaling types.
scalingType
Output format types.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const word canonicalName
The canonical name ("lumpedPointMovement") for the dictionary.
Templated 2D Barycentric derived from VectorSpace. Has 3 components, one of which is redundant...
Definition: Barycentric2D.H:50
static const Enum< inputFormatType > formatNames
Names for the input format types.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:342
A simple linear interpolator between two locations, which are referenced by index.
dynamicFvMesh & mesh
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:113
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:65
bool forcesAndMoments(const polyMesh &pmesh, List< vector > &forces, List< vector > &moments) const
The forces and moments acting on each pressure-zone.
void clear()
Remove all entries from table.
Definition: HashTable.C:725
bool hasInterpolator(const pointPatch &fpatch) const
Check if patch control exists for specified patch.
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
outputFormatType
Output format types.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual const word & name() const =0
Return name.
End list [isseparator].
Definition: token.H:162
Vector< scalar > vector
Definition: vector.H:57
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
label size() const noexcept
The number of elements in the list.
const dimensionSet dimPressure
bool test(const label pos) const
Test for True value at specified position, never auto-vivify entries.
Definition: bitSet.H:326
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTable.H:337
void writeList(vtk::formatter &fmt, const UList< uint8_t > &values)
Write a list of uint8_t values.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging 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...
tmp< pointField > pointsDisplacement(const pointPatch &fpatch, const pointField &points0) const
Displace points according to the current state.
OBJstream os(runTime.globalPath()/outputName)
bool readState()
Read state from file, applying relaxation as requested.
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:405
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Find and return a sub-dictionary as a copy, otherwise return an empty dictionary. ...
Definition: dictionary.C:521
void scalePoints(const scalar scaleFactor)
Scale points by given factor.
void writeDict(Ostream &os) const
Write axis, locations, division as a dictionary.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
Point vec() const
Return start-to-end vector.
Definition: lineI.H:109
void checkPatchControl(const polyPatch &pp) const
Check if patch control exists for specified patch.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))
bool couplingPending(const label timeIndex) const
Check if coupling is pending (according to the calcFrequency)
Non-pointer based hierarchical recursive searching.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
#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
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
virtual label index() const =0
Return the index of this patch in the pointBoundaryMesh.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
const tensorField & rotations() const
The local-to-global transformation for each point.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:139
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
const T2 & second() const noexcept
Access the second element.
Definition: Tuple2.H:142
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
const polyBoundaryMesh & patches
const pointField & points() const
The points corresponding to mass centres.
static int debug
Debug switch.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
List< label > sortedToc(const UList< bool > &bools)
Return the (sorted) values corresponding to &#39;true&#39; entries.
Definition: BitOps.C:195
messageStream Info
Information stream (stdout output on master, null elsewhere)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:124
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...
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual Ostream & endEntry()
Write end entry (&#39;;&#39;) followed by newline.
Definition: Ostream.C:117
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
const T & second() const noexcept
Access the second element.
Definition: Pair.H:147
const T1 & first() const noexcept
Access the first element.
Definition: Tuple2.H:132
bool readData(Istream &is, const quaternion::eulerOrder rotOrder=quaternion::eulerOrder::ZXZ, const bool degrees=false)
Read input as dictionary content.
Tensor of scalars, i.e. Tensor<scalar>.
bool hasPatchControl(const label patchIndex) const
Check if patch control exists for specified patch.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
lumpedPointMovement()
Default construct.
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
The state of lumped points corresponds to positions and rotations.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
virtual const labelList & meshPoints() const =0
Return mesh points.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
label timeIndex
Definition: getTimeIndex.H:24
const dimensionSet & dimensions() const noexcept
Return dimensions.
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.
A topoSetPointSource to select all points based on usage in given faceSet(s).
Definition: faceToPoint.H:170
virtual label size() const =0
Return size.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
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