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