Loading...
Searching...
No Matches
particle.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) 2011-2017, 2020 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "particle.H"
30#include "transform.H"
31#include "treeDataCell.H"
32#include "cubicEqn.H"
34#include "indexedOctree.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
41}
42
43const Foam::label Foam::particle::maxNBehind_ = 10;
44
45Foam::label Foam::particle::particleCount_ = 0;
46
48
50(
51 Foam::debug::infoSwitch("writeLagrangianPositions", 1)
52);
53
55(
56 "writeLagrangianPositions",
57 bool,
59);
60
61
62// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
63
64void Foam::particle::stationaryTetReverseTransform
65(
66 vector& centre,
67 scalar& detA,
68 barycentricTensor& T
69) const
70{
71 barycentricTensor A = stationaryTetTransform();
72
73 const vector ab = A.b() - A.a();
74 const vector ac = A.c() - A.a();
75 const vector ad = A.d() - A.a();
76 const vector bc = A.c() - A.b();
77 const vector bd = A.d() - A.b();
78
79 centre = A.a();
80
81 detA = ab & (ac ^ ad);
82
84 (
85 bd ^ bc,
86 ac ^ ad,
87 ad ^ ab,
88 ab ^ ac
89 );
90}
91
92
93void Foam::particle::movingTetReverseTransform
94(
95 const scalar fraction,
96 Pair<vector>& centre,
99) const
100{
101 Pair<barycentricTensor> A = movingTetTransform(fraction);
102
103 const Pair<vector> ab(A[0].b() - A[0].a(), A[1].b() - A[1].a());
104 const Pair<vector> ac(A[0].c() - A[0].a(), A[1].c() - A[1].a());
105 const Pair<vector> ad(A[0].d() - A[0].a(), A[1].d() - A[1].a());
106 const Pair<vector> bc(A[0].c() - A[0].b(), A[1].c() - A[1].b());
107 const Pair<vector> bd(A[0].d() - A[0].b(), A[1].d() - A[1].b());
108
109 centre[0] = A[0].a();
110 centre[1] = A[1].a();
111
112 detA[0] = ab[0] & (ac[0] ^ ad[0]);
113 detA[1] =
114 (ab[1] & (ac[0] ^ ad[0]))
115 + (ab[0] & (ac[1] ^ ad[0]))
116 + (ab[0] & (ac[0] ^ ad[1]));
117 detA[2] =
118 (ab[0] & (ac[1] ^ ad[1]))
119 + (ab[1] & (ac[0] ^ ad[1]))
120 + (ab[1] & (ac[1] ^ ad[0]));
121 detA[3] = ab[1] & (ac[1] ^ ad[1]);
122
124 (
125 bd[0] ^ bc[0],
126 ac[0] ^ ad[0],
127 ad[0] ^ ab[0],
128 ab[0] ^ ac[0]
129 );
131 (
132 (bd[0] ^ bc[1]) + (bd[1] ^ bc[0]),
133 (ac[0] ^ ad[1]) + (ac[1] ^ ad[0]),
134 (ad[0] ^ ab[1]) + (ad[1] ^ ab[0]),
135 (ab[0] ^ ac[1]) + (ab[1] ^ ac[0])
136 );
138 (
139 bd[1] ^ bc[1],
140 ac[1] ^ ad[1],
141 ad[1] ^ ab[1],
142 ab[1] ^ ac[1]
143 );
144}
145
146
147void Foam::particle::reflect()
148{
149 std::swap(coordinates_.c(), coordinates_.d());
150}
151
152
153void Foam::particle::rotate(const bool reverse)
154{
155 if (!reverse)
156 {
157 scalar temp = coordinates_.b();
158 coordinates_.b() = coordinates_.c();
159 coordinates_.c() = coordinates_.d();
160 coordinates_.d() = temp;
161 }
162 else
163 {
164 scalar temp = coordinates_.d();
165 coordinates_.d() = coordinates_.c();
166 coordinates_.c() = coordinates_.b();
167 coordinates_.b() = temp;
168 }
169}
170
171
172void Foam::particle::changeTet(const label tetTriI)
173{
174 const bool isOwner = mesh_.faceOwner()[tetFacei_] == celli_;
175
176 const label firstTetPtI = 1;
177 const label lastTetPtI = mesh_.faces()[tetFacei_].size() - 2;
178
179 if (tetTriI == 1)
180 {
181 changeFace(tetTriI);
182 }
183 else if (tetTriI == 2)
184 {
185 if (isOwner)
186 {
187 if (tetPti_ == lastTetPtI)
188 {
189 changeFace(tetTriI);
190 }
191 else
192 {
193 reflect();
194 tetPti_ += 1;
195 }
196 }
197 else
198 {
199 if (tetPti_ == firstTetPtI)
200 {
201 changeFace(tetTriI);
202 }
203 else
204 {
205 reflect();
206 tetPti_ -= 1;
207 }
208 }
209 }
210 else if (tetTriI == 3)
211 {
212 if (isOwner)
213 {
214 if (tetPti_ == firstTetPtI)
215 {
216 changeFace(tetTriI);
217 }
218 else
219 {
220 reflect();
221 tetPti_ -= 1;
222 }
223 }
224 else
225 {
226 if (tetPti_ == lastTetPtI)
227 {
228 changeFace(tetTriI);
229 }
230 else
231 {
232 reflect();
233 tetPti_ += 1;
234 }
235 }
236 }
237 else
238 {
240 << "Changing tet without changing cell should only happen when the "
241 << "track is on triangle 1, 2 or 3."
242 << exit(FatalError);
243 }
244}
245
246
247void Foam::particle::changeFace(const label tetTriI)
248{
249 // Get the old topology
250 const triFace triOldIs(currentTetIndices().faceTriIs(mesh_));
251
252 // Get the shared edge and the pre-rotation
253 edge sharedEdge;
254 if (tetTriI == 1)
255 {
256 sharedEdge = edge(triOldIs[1], triOldIs[2]);
257 }
258 else if (tetTriI == 2)
259 {
260 sharedEdge = edge(triOldIs[2], triOldIs[0]);
261 }
262 else if (tetTriI == 3)
263 {
264 sharedEdge = edge(triOldIs[0], triOldIs[1]);
265 }
266 else
267 {
269 << "Changing face without changing cell should only happen when the"
270 << " track is on triangle 1, 2 or 3."
271 << exit(FatalError);
272
273 sharedEdge = edge(-1, -1);
274 }
275
276 // Find the face in the same cell that shares the edge, and the
277 // corresponding tetrahedra point
278 tetPti_ = -1;
279 for (const label newFaceI : mesh_.cells()[celli_])
280 {
281 const Foam::face& newFace = mesh_.faces()[newFaceI];
282 const label newOwner = mesh_.faceOwner()[newFaceI];
283
284 // Exclude the current face
285 if (tetFacei_ == newFaceI)
286 {
287 continue;
288 }
289
290 // Loop over the edges, looking for the shared one. Note that we have to
291 // match the direction of the edge as well as the end points in order to
292 // avoid false positives when dealing with coincident ACMI faces.
293 const label edgeComp = newOwner == celli_ ? -1 : +1;
294 label edgeI = 0;
295 for
296 (
297 ;
298 edgeI < newFace.size()
299 && edge::compare(sharedEdge, newFace.edge(edgeI)) != edgeComp;
300 ++ edgeI
301 );
302
303 // If the face does not contain the edge, then move on to the next face
304 if (edgeI >= newFace.size())
305 {
306 continue;
307 }
308
309 // Make the edge index relative to the base point
310 const label newBaseI = max(0, mesh_.tetBasePtIs()[newFaceI]);
311 edgeI = (edgeI - newBaseI + newFace.size()) % newFace.size();
312
313 // If the edge is next the base point (i.e., the index is 0 or n - 1),
314 // then we swap it for the adjacent edge. This new edge is opposite the
315 // base point, and defines the tet with the original edge in it.
316 edgeI = min(max(1, edgeI), newFace.size() - 2);
317
318 // Set the new face and tet point
319 tetFacei_ = newFaceI;
320 tetPti_ = edgeI;
321
322 // Exit the loop now that the tet point has been found
323 break;
324 }
325
326 if (tetPti_ == -1)
327 {
329 << "The search for an edge-connected face and tet-point failed."
330 << exit(FatalError);
331 }
332
333 // Pre-rotation puts the shared edge opposite the base of the tetrahedron
334 if (sharedEdge.otherVertex(triOldIs[1]) == -1)
335 {
336 rotate(false);
337 }
338 else if (sharedEdge.otherVertex(triOldIs[2]) == -1)
339 {
340 rotate(true);
341 }
342
343 // Get the new topology
344 const triFace triNewIs = currentTetIndices().faceTriIs(mesh_);
345
346 // Reflect to account for the change of triangle orientation on the new face
347 reflect();
348
349 // Post rotation puts the shared edge back in the correct location
350 if (sharedEdge.otherVertex(triNewIs[1]) == -1)
351 {
352 rotate(true);
353 }
354 else if (sharedEdge.otherVertex(triNewIs[2]) == -1)
355 {
356 rotate(false);
357 }
358}
359
360
361void Foam::particle::changeCell()
362{
363 // Set the cell to be the one on the other side of the face
364 const label ownerCellI = mesh_.faceOwner()[tetFacei_];
365 const bool isOwner = celli_ == ownerCellI;
366 celli_ = isOwner ? mesh_.faceNeighbour()[tetFacei_] : ownerCellI;
367 // Reflect to account for the change of triangle orientation in the new cell
368 reflect();
369}
370
371
372void Foam::particle::changeToMasterPatch()
373{
374 label thisPatch = patch();
375
376 for (const label otherFaceI : mesh_.cells()[celli_])
377 {
378 // Skip the current face and any internal faces
379 if (facei_ == otherFaceI || mesh_.isInternalFace(otherFaceI))
380 {
381 continue;
382 }
383
384 // Compare the two faces. If they are the same, chose the one with the
385 // lower patch index. In the case of an ACMI-wall pair, this will be
386 // the ACMI, which is what we want.
387 const Foam::face& thisFace = mesh_.faces()[facei_];
388 const Foam::face& otherFace = mesh_.faces()[otherFaceI];
389 if (face::compare(thisFace, otherFace) != 0)
390 {
391 const label otherPatch =
392 mesh_.boundaryMesh().whichPatch(otherFaceI);
393 if (thisPatch > otherPatch)
394 {
395 facei_ = otherFaceI;
396 thisPatch = otherPatch;
397 }
398 }
399 }
400
401 tetFacei_ = facei_;
402}
403
404
405void Foam::particle::locate
406(
407 const vector& position,
408 const vector* direction,
409 const label celli,
410 const bool boundaryFail,
411 const string& boundaryMsg
412)
413{
414 if (debug)
415 {
416 Pout<< "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
417 }
418
419 celli_ = celli;
420
421 // Find the cell, if it has not been given
422 if (celli_ < 0)
423 {
424 celli_ = mesh_.cellTree().findInside(position);
425 if (celli_ < 0)
426 {
428 << "Cell not found for particle position " << position << "."
429 << exit(FatalError);
430 }
431 }
432
433 // Perturbing displacement to avoid zero in case position is the
434 // cellCentre
435 const vector displacement =
436 position
437 - mesh_.cellCentres()[celli_]
438 + vector(VSMALL, VSMALL, VSMALL);
439
440 // Loop all cell tets to find the one containing the position. Track
441 // through each tet from the cell centre. If a tet contains the position
442 // then the track will end with a single trackToTri.
443 const Foam::cell& c = mesh_.cells()[celli_];
444 scalar minF = VGREAT;
445 label minTetFacei = -1, minTetPti = -1;
446 forAll(c, cellTetFacei)
447 {
448 const Foam::face& f = mesh_.faces()[c[cellTetFacei]];
449 for (label tetPti = 1; tetPti < f.size() - 1; ++tetPti)
450 {
451 coordinates_ = barycentric(1, 0, 0, 0);
452 tetFacei_ = c[cellTetFacei];
453 tetPti_ = tetPti;
454 facei_ = -1;
455
456 label tetTriI = -1;
457 const scalar f = trackToTri(displacement, 0, tetTriI);
458
459 if (tetTriI == -1)
460 {
461 return;
462 }
463
464 if (f < minF)
465 {
466 minF = f;
467 minTetFacei = tetFacei_;
468 minTetPti = tetPti_;
469 }
470 }
471 }
472
473 // The particle must be (hopefully only slightly) outside the cell. Track
474 // into the tet which got the furthest.
475 coordinates_ = barycentric(1, 0, 0, 0);
476 tetFacei_ = minTetFacei;
477 tetPti_ = minTetPti;
478 facei_ = -1;
479 track(displacement, 0);
480 if (!onFace())
481 {
482 return;
483 }
484
485 // If we are here then we hit a boundary
486 if (boundaryFail)
487 {
488 FatalErrorInFunction << boundaryMsg << exit(FatalError);
489 }
490 else
491 {
492 static label nWarnings = 0;
493 static const label maxNWarnings = 100;
494 if ((nWarnings < maxNWarnings) && boundaryFail)
495 {
496 WarningInFunction << boundaryMsg.c_str() << endl;
497 ++ nWarnings;
498 }
499 if (nWarnings == maxNWarnings)
500 {
502 << "Suppressing any further warnings about particles being "
503 << "located outside of the mesh." << endl;
504 ++ nWarnings;
505 }
507
508}
509
510
511// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
512
514(
515 const polyMesh& mesh,
517 const label celli,
518 const label tetFacei,
519 const label tetPti
520)
521:
522 mesh_(mesh),
523 coordinates_(coordinates),
524 celli_(celli),
525 tetFacei_(tetFacei),
526 tetPti_(tetPti),
527 facei_(-1),
528 stepFraction_(1.0),
529 behind_(0.0),
530 nBehind_(0),
531 origProc_(Pstream::myProcNo()),
532 origId_(getNewParticleID())
533{}
534
535
537(
538 const polyMesh& mesh,
539 const vector& position,
540 const label celli
541)
542:
543 mesh_(mesh),
544 coordinates_(-VGREAT, -VGREAT, -VGREAT, -VGREAT),
545 celli_(celli),
546 tetFacei_(-1),
547 tetPti_(-1),
548 facei_(-1),
549 stepFraction_(1.0),
550 behind_(0.0),
551 nBehind_(0),
552 origProc_(Pstream::myProcNo()),
553 origId_(getNewParticleID())
554{
555 locate
556 (
557 position,
558 nullptr,
559 celli,
560 false,
561 "Particle initialised with a location outside of the mesh"
562 );
563}
564
565
567(
568 const polyMesh& mesh,
569 const vector& position,
570 const label celli,
571 const label tetFacei,
572 const label tetPti,
573 bool doLocate
574)
575:
576 mesh_(mesh),
577 coordinates_(-VGREAT, -VGREAT, -VGREAT, -VGREAT),
578 celli_(celli),
579 tetFacei_(tetFacei),
580 tetPti_(tetPti),
581 facei_(-1),
582 stepFraction_(1.0),
583 behind_(0.0),
584 nBehind_(0),
585 origProc_(Pstream::myProcNo()),
586 origId_(getNewParticleID())
587{
588 if (doLocate)
589 {
590 locate
591 (
592 position,
593 nullptr,
594 celli,
595 false,
596 "Particle initialised with a location outside of the mesh"
597 );
598 }
599}
600
601
602Foam::particle::particle(const particle& p, const polyMesh& mesh)
603:
604 mesh_(mesh),
605 coordinates_(p.coordinates_),
606 celli_(p.celli_),
607 tetFacei_(p.tetFacei_),
608 tetPti_(p.tetPti_),
609 facei_(p.facei_),
610 stepFraction_(p.stepFraction_),
611 behind_(p.behind_),
612 nBehind_(p.nBehind_),
613 origProc_(p.origProc_),
614 origId_(p.origId_)
615{}
616
617
620 particle(p, p.mesh())
621{}
622
623
624// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
625
626Foam::scalar Foam::particle::track
627(
628 const vector& displacement,
629 const scalar fraction
630)
631{
632 scalar f = trackToFace(displacement, fraction);
633
634 while (onInternalFace())
635 {
636 changeCell();
637
638 f *= trackToFace(f*displacement, f*fraction);
639 }
640
641 return f;
642}
643
644
646(
647 const vector& displacement,
648 const scalar fraction
649)
650{
651 scalar f = 1;
652
653 label tetTriI = onFace() ? 0 : -1;
654
655 facei_ = -1;
656
657 // Loop the tets in the current cell
658 while (nBehind_ < maxNBehind_)
659 {
660 f *= trackToTri(f*displacement, f*fraction, tetTriI);
661
662 if (tetTriI == -1)
663 {
664 // The track has completed within the current tet
665 return 0;
666 }
667 else if (tetTriI == 0)
668 {
669 // The track has hit a face, so set the current face and return
670 facei_ = tetFacei_;
671 return f;
672 }
673 else
674 {
675 // Move to the next tet and continue the track
676 changeTet(tetTriI);
677 }
678 }
679
680 // Warn if stuck, and incorrectly advance the step fraction to completion
681 static label stuckID = -1, stuckProc = -1;
682 if (origId_ != stuckID && origProc_ != stuckProc)
683 {
685 << "Particle #" << origId_ << " got stuck at " << position()
686 << endl;
687 }
688
689 stuckID = origId_;
690 stuckProc = origProc_;
691
692 stepFraction_ += f*fraction;
693
694 behind_ = 0;
695 nBehind_ = 0;
696
697 return 0;
698}
699
700
702(
703 const vector& displacement,
704 const scalar fraction,
705 label& tetTriI
706)
707{
708 const vector x0 = position();
709 const vector x1 = displacement;
710 const barycentric y0 = coordinates_;
711
712 if (debug)
713 {
714 Pout<< "Particle " << origId() << endl << "Tracking from " << x0
715 << " along " << x1 << " to " << x0 + x1 << endl;
716 }
717
718 // Get the tet geometry
719 vector centre;
720 scalar detA;
722 stationaryTetReverseTransform(centre, detA, T);
723
724 if (debug)
725 {
726 Pout<< "Tet points: " << currentTetIndices().tet(mesh_) << endl
727 << "Tet determinant = " << detA << endl
728 << "Start local coordinates = " << y0 << endl;
729 }
730
731 // Calculate the local tracking displacement
732 barycentric Tx1(x1 & T);
733
734 if (debug)
735 {
736 Pout<< "Local displacement = " << Tx1 << "/" << detA << endl;
737 }
738
739 // Calculate the hit fraction
740 label iH = -1;
741 scalar muH = detA <= 0 ? VGREAT : 1/detA;
742 for (label i = 0; i < 4; ++ i)
743 {
744 if (Tx1[i] < - detA*SMALL)
745 {
746 scalar mu = - y0[i]/Tx1[i];
747
748 if (debug)
749 {
750 Pout<< "Hit on tet face " << i << " at local coordinate "
751 << y0 + mu*Tx1 << ", " << mu*detA*100 << "% of the "
752 << "way along the track" << endl;
753 }
754
755 if (0 <= mu && mu < muH)
756 {
757 iH = i;
758 muH = mu;
759 }
760 }
761 }
762
763 // Set the new coordinates
764 barycentric yH = y0 + muH*Tx1;
765
766 // Clamp to zero any negative coordinates generated by round-off error
767 for (label i = 0; i < 4; ++ i)
768 {
769 yH.replace(i, i == iH ? 0 : max(0, yH[i]));
770 }
771
772 // Re-normalise if within the tet
773 if (iH == -1)
774 {
775 yH /= cmptSum(yH);
776 }
777
778 // Set the new position and hit index
779 coordinates_ = yH;
780 tetTriI = iH;
781
782 if (debug)
783 {
784 if (iH != -1)
785 {
786 Pout<< "Track hit tet face " << iH << " first" << endl;
787 }
788 else
789 {
790 Pout<< "Track hit no tet faces" << endl;
791 }
792 Pout<< "End local coordinates = " << yH << endl
793 << "End global coordinates = " << position() << endl
794 << "Tracking displacement = " << position() - x0 << endl
795 << muH*detA*100 << "% of the step from " << stepFraction_ << " to "
796 << stepFraction_ + fraction << " completed" << endl << endl;
797 }
798
799 // Set the proportion of the track that has been completed
800 stepFraction_ += fraction*muH*detA;
801
802 if (debug)
803 {
804 Pout << "Step Fraction : " << stepFraction_ << endl;
805 Pout << "fraction*muH*detA : " << fraction*muH*detA << endl;
806 Pout << "static muH : " << muH << endl;
807 Pout << "origId() : " << origId() << endl;
808 }
809
810 // Accumulate displacement behind
811 if (detA <= 0 || nBehind_ > 0)
812 {
813 behind_ += muH*detA*mag(displacement);
814
815 if (behind_ > 0)
816 {
817 behind_ = 0;
818 nBehind_ = 0;
819 }
820 else
821 {
822 ++ nBehind_;
824 }
825
826 return iH != -1 ? 1 - muH*detA : 0;
827}
828
829
831(
832 const vector& displacement,
833 const scalar fraction,
834 label& tetTriI
835)
836{
837 const vector x0 = position();
838 const vector x1 = displacement;
839 const barycentric y0 = coordinates_;
840
841 if (debug)
842 {
843 Pout<< "Particle " << origId() << endl << "Tracking from " << x0
844 << " along " << x1 << " to " << x0 + x1 << endl;
845 }
846
847 // Get the tet geometry
848 Pair<vector> centre;
851 movingTetReverseTransform(fraction, centre, detA, T);
852
853 if (debug)
854 {
855 Pair<vector> o, b, v1, v2;
856 movingTetGeometry(fraction, o, b, v1, v2);
857 Pout<< "Tet points o=" << o[0] << ", b=" << b[0]
858 << ", v1=" << v1[0] << ", v2=" << v2[0] << endl
859 << "Tet determinant = " << detA[0] << endl
860 << "Start local coordinates = " << y0[0] << endl;
861 }
862
863
864 // Get the relative global position
865 const vector x0Rel = x0 - centre[0];
866 const vector x1Rel = x1 - centre[1];
867
868 // Form the determinant and hit equations
869 cubicEqn detAEqn(sqr(detA[0])*detA[3], detA[0]*detA[2], detA[1], 1);
870 const barycentric yC(1, 0, 0, 0);
871 const barycentric hitEqnA =
872 ((x1Rel & T[2]) + detA[3]*yC)*sqr(detA[0]);
873 const barycentric hitEqnB =
874 ((x1Rel & T[1]) + (x0Rel & T[2]) + detA[2]*yC)*detA[0];
875 const barycentric hitEqnC =
876 ((x1Rel & T[0]) + (x0Rel & T[1]) + detA[1]*yC);
877 const barycentric hitEqnD = y0;
879 forAll(hitEqn, i)
880 {
881 hitEqn[i] = cubicEqn(hitEqnA[i], hitEqnB[i], hitEqnC[i], hitEqnD[i]);
882 }
883
884 if (debug)
885 {
886 for (label i = 0; i < 4; ++ i)
887 {
888 Pout<< (i ? " " : "Hit equation ") << i << " = "
889 << hitEqn[i] << endl;
890 }
891 Pout<< " DetA equation = " << detA << endl;
892 }
893
894 // Calculate the hit fraction
895 label iH = -1;
896 scalar muH = detA[0] <= 0 ? VGREAT : 1/detA[0];
897 for (label i = 0; i < 4; ++i)
898 {
899 const Roots<3> mu = hitEqn[i].roots();
900
901 for (label j = 0; j < 3; ++j)
902 {
903 if
904 (
905 mu.type(j) == roots::real
906 && hitEqn[i].derivative(mu[j]) < - detA[0]*SMALL
907 )
908 {
909 if (debug)
910 {
911 const barycentric yH
912 (
913 hitEqn[0].value(mu[j]),
914 hitEqn[1].value(mu[j]),
915 hitEqn[2].value(mu[j]),
916 hitEqn[3].value(mu[j])
917 );
918 const scalar detAH = detAEqn.value(mu[j]);
919
920 Pout<< "Hit on tet face " << i << " at local coordinate "
921 << (std::isnormal(detAH) ? name(yH/detAH) : "???")
922 << ", " << mu[j]*detA[0]*100 << "% of the "
923 << "way along the track" << endl;
924
925 Pout<< "derivative : " << hitEqn[i].derivative(mu[j]) << nl
926 << " coord " << j << " mu[j]: " << mu[j] << nl
927 << " hitEq " << i << endl;
928 }
929
930 if (0 <= mu[j] && mu[j] < muH)
931 {
932 iH = i;
933 muH = mu[j];
934 }
935 }
936 }
937 }
938
939 // Set the new coordinates
940 barycentric yH
941 (
942 hitEqn[0].value(muH),
943 hitEqn[1].value(muH),
944 hitEqn[2].value(muH),
945 hitEqn[3].value(muH)
946 );
947 // !!! <-- This fails if the tet collapses onto the particle, as detA tends
948 // to zero at the hit. In this instance, we can differentiate the hit and
949 // detA polynomials to find a limiting location, but this will not be on a
950 // triangle. We will then need to do a second track through the degenerate
951 // tet to find the final hit position. This second track is over zero
952 // distance and therefore can be of the static mesh type. This has not yet
953 // been implemented.
954 const scalar detAH = detAEqn.value(muH);
955 if (!std::isnormal(detAH))
956 {
958 << "A moving tet collapsed onto a particle. This is not supported. "
959 << "The mesh is too poor, or the motion too severe, for particle "
960 << "tracking to function." << exit(FatalError);
961 }
962 yH /= detAH;
963
964 // Clamp to zero any negative coordinates generated by round-off error
965 for (label i = 0; i < 4; ++ i)
966 {
967 yH.replace(i, i == iH ? 0 : max(0, yH[i]));
968 }
969
970 // Re-normalise if within the tet
971 if (iH == -1)
972 {
973 yH /= cmptSum(yH);
974 }
975
976 // Set the new position and hit index
977 coordinates_ = yH;
978 tetTriI = iH;
979
980 scalar advance = muH*detA[0];
981
982 // Set the proportion of the track that has been completed
983 stepFraction_ += fraction*advance;
984
985 // Accumulate displacement behind
986 if (detA[0] <= 0 || nBehind_ > 0)
987 {
988 behind_ += muH*detA[0]*mag(displacement);
989
990 if (behind_ > 0)
991 {
992 behind_ = 0;
993 nBehind_ = 0;
994 }
995 else
996 {
997 ++ nBehind_;
998 }
999 }
1000
1001 if (debug)
1002 {
1003 if (iH != -1)
1004 {
1005 Pout<< "Track hit tet face " << iH << " first" << endl;
1006 }
1007 else
1008 {
1009 Pout<< "Track hit no tet faces" << endl;
1010 }
1011// Pout<< "End local coordinates = " << yH << endl
1012// << "End global coordinates = " << position() << endl
1013// << "Tracking displacement = " << position() - x0 << endl
1014// << muH*detA[0]*100 << "% of the step from " << stepFraction_
1015// << " to " << stepFraction_ + fraction << " completed" << endl
1016// << endl;
1018
1019
1020 return iH != -1 ? 1 - muH*detA[0] : 0;
1021}
1022
1023
1024Foam::scalar Foam::particle::trackToTri
1025(
1026 const vector& displacement,
1027 const scalar fraction,
1028 label& tetTriI
1029)
1030{
1031 if ((mesh_.moving() && (stepFraction_ != 1 || fraction != 0)))
1032 {
1033 return trackToMovingTri(displacement, fraction, tetTriI);
1034 }
1035 else
1036 {
1037 return trackToStationaryTri(displacement, fraction, tetTriI);
1038 }
1039}
1040
1041
1043{
1044 if (cmptMin(mesh_.geometricD()) == -1)
1045 {
1046 vector pos = position(), posC = pos;
1048 return pos - posC;
1049 }
1050 else
1051 {
1052 return vector::zero;
1053 }
1055
1056
1059
1060
1062{}
1063
1064
1066{
1067 // Convert the face index to be local to the processor patch
1068 facei_ = mesh_.boundaryMesh()[patch()].whichFace(facei_);
1069}
1070
1071
1073(
1074 const label patchi,
1075 trackingData& td
1076)
1077{
1078 const coupledPolyPatch& ppp =
1079 refCast<const coupledPolyPatch>(mesh_.boundaryMesh()[patchi]);
1080
1081 if (!ppp.parallel())
1082 {
1083 const tensor& T =
1084 (
1085 ppp.forwardT().size() == 1
1086 ? ppp.forwardT()[0]
1087 : ppp.forwardT()[facei_]
1088 );
1089 transformProperties(T);
1090 }
1091 else if (ppp.separated())
1092 {
1093 const vector& s =
1094 (
1095 (ppp.separation().size() == 1)
1096 ? ppp.separation()[0]
1097 : ppp.separation()[facei_]
1098 );
1099 transformProperties(-s);
1100 }
1101
1102 // Set the topology
1103 celli_ = ppp.faceCells()[facei_];
1104 facei_ += ppp.start();
1105 tetFacei_ = facei_;
1106 // Faces either side of a coupled patch are numbered in opposite directions
1107 // as their normals both point away from their connected cells. The tet
1108 // point therefore counts in the opposite direction from the base point.
1109 tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;
1110
1111 // Reflect to account for the change of triangle orientation in the new cell
1112 reflect();
1113
1114 // Note that the position does not need transforming explicitly. The face-
1115 // triangle on the receive patch is the transformation of the one on the
1116 // send patch, so whilst the barycentric coordinates remain the same, the
1117 // change of triangle implicitly transforms the position.
1118}
1119
1120
1122(
1124)
1125{
1126 // Get the transformed position
1127 const vector pos = transform.invTransformPosition(position());
1128
1129 // Break the topology
1130 celli_ = -1;
1131 tetFacei_ = -1;
1132 tetPti_ = -1;
1133 facei_ = -1;
1134
1135 // Store the position in the barycentric data
1136 coordinates_ = barycentric(1 - cmptSum(pos), pos.x(), pos.y(), pos.z());
1137
1138 // Transform the properties
1139 transformProperties(- transform.t());
1140 if (transform.hasR())
1141 {
1143 }
1144}
1145
1146
1148{
1149 // Get the position from the barycentric data
1150 const vector pos(coordinates_.b(), coordinates_.c(), coordinates_.d());
1151
1152 // Create some arbitrary topology for the supplied cell
1153 celli_ = celli;
1154 tetFacei_ = mesh_.cells()[celli_][0];
1155 tetPti_ = 1;
1156 facei_ = -1;
1157
1158 // Get the reverse transform and directly set the coordinates from the
1159 // position. This isn't likely to be correct; the particle is probably not
1160 // in this tet. It will, however, generate the correct vector when the
1161 // position method is called. A referred particle should never be tracked,
1162 // so this approximate topology is good enough. By using the nearby cell we
1163 // minimise the error associated with the incorrect topology.
1164 coordinates_ = barycentric(1, 0, 0, 0);
1165 if (mesh_.moving() && stepFraction_ != 1)
1166 {
1167 Pair<vector> centre;
1168 FixedList<scalar, 4> detA;
1169 FixedList<barycentricTensor, 3> T;
1170 movingTetReverseTransform(0, centre, detA, T);
1171 coordinates_ += (pos - centre[0]) & T[0]/detA[0];
1172 }
1173 else
1174 {
1175 vector centre;
1176 scalar detA;
1178 stationaryTetReverseTransform(centre, detA, T);
1179 coordinates_ += (pos - centre) & T/detA;
1180 }
1181}
1182
1183
1184Foam::label Foam::particle::procTetPt
1185(
1186 const polyMesh& procMesh,
1187 const label procCell,
1188 const label procTetFace
1189) const
1190{
1191 // The tet point on the procMesh differs from the current tet point if the
1192 // mesh and procMesh faces are of differing orientation. The change is the
1193 // same as in particle::correctAfterParallelTransfer.
1194
1195 if
1196 (
1197 (mesh_.faceOwner()[tetFacei_] == celli_)
1198 == (procMesh.faceOwner()[procTetFace] == procCell)
1199 )
1200 {
1201 return tetPti_;
1202 }
1203 else
1204 {
1205 return procMesh.faces()[procTetFace].size() - 1 - tetPti_;
1206 }
1207}
1208
1209
1211(
1212 const vector& position,
1213 const mapPolyMesh& mapper
1214)
1215{
1216 locate
1217 (
1218 position,
1219 nullptr,
1220 mapper.reverseCellMap()[celli_],
1221 true,
1222 "Particle mapped to a location outside of the mesh"
1223 );
1224}
1225
1226
1227void Foam::particle::relocate(const point& position, const label celli)
1228{
1229 locate
1230 (
1231 position,
1232 nullptr,
1233 celli,
1234 true,
1235 "Particle mapped to a location outside of the mesh"
1236 );
1237}
1238
1239
1240// * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
1242bool Foam::operator==(const particle& pA, const particle& pB)
1243{
1244 return (pA.origProc() == pB.origProc() && pA.origId() == pB.origId());
1245}
1246
1247
1248bool Foam::operator!=(const particle& pA, const particle& pB)
1249{
1250 return !(pA == pB);
1251}
1252
1253
1254// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Inter-processor communications stream.
Definition Pstream.H:59
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition Roots.H:71
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
void replace(const direction, const Cmpt &)
The coupledPolyPatch is an abstract base class for patches that couple regions of the computational d...
virtual bool separated() const
Are the planes separated.
virtual bool parallel() const
Are the cyclic planes parallel.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Container to encapsulate various operations for cubic equation of the forms with real coefficients:
Definition cubicEqn.H:111
scalar value(const scalar x) const
Evaluate the cubic equation at x.
Definition cubicEqnI.H:46
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
static int compare(const edge &a, const edge &b)
Compare edges.
Definition edgeI.H:26
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition faceI.H:158
static int compare(const face &a, const face &b)
Compare faces.
Definition face.C:273
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseCellMap() const noexcept
Reverse cell map.
Base particle class.
Definition particle.H:72
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
Definition particle.C:1140
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition particle.C:1050
label procTetPt(const polyMesh &procMesh, const label procCell, const label procTetFace) const
Return the tet point appropriate for decomposition or reconstruction.
Definition particle.C:1178
tetIndices currentTetIndices() const noexcept
Return indices of the current tet that the particle occupies.
Definition particleI.H:236
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition particle.C:1035
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition particle.C:695
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition particle.C:1066
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but also stops on internal faces.
Definition particle.C:639
vector position() const
Return current particle position.
Definition particleI.H:283
void prepareForInteractionListReferral(const vectorTensorTransform &transform)
Break the topology and store the particle position so that the.
Definition particle.C:1115
bool onFace() const noexcept
Is the particle on a face?
Definition particleI.H:259
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition particle.C:1204
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition particleI.H:116
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition particle.C:824
label patch() const
Return the index of patch that the particle is on.
Definition particleI.H:277
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
Definition particle.C:1058
static label particleCount_
Cumulative particle counter - used to provide unique ID.
Definition particle.H:460
static bool writeLagrangianPositions
Write particle positions file (v1706 format and earlier) Default is true (disable in etc/controlDict)...
Definition particle.H:472
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition particle.C:620
bool onInternalFace() const noexcept
Is the particle on an internal face?
Definition particleI.H:265
label getNewParticleID() const
Get unique particle creation id.
Definition particleI.H:96
label origId() const noexcept
Return the particle ID on the originating processor.
Definition particleI.H:194
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
void relocate(const point &position, const label celli=-1)
Set the addressing based on the provided position.
Definition particle.C:1220
scalar trackToTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but also stops on tet triangles. On.
Definition particle.C:1018
static bool writeLagrangianCoordinates
Write particle coordinates file (v1712 and later) Default is true.
Definition particle.H:466
label origProc() const noexcept
Return the originating processor ID.
Definition particleI.H:182
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
const labelUList & faceCells() const
Return face-cell addressing.
Definition polyPatch.C:401
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition triFace.H:68
Vector-tensor class used to perform translations and rotations in 3D space.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
int infoSwitch(const char *name, const int deflt=0)
Lookup info switch or add default value.
Definition debug.C:228
const std::string patch
OpenFOAM patch number as a std::string.
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition meshTools.C:622
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition meshTools.C:548
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
dimensionedScalar pos(const dimensionedScalar &ds)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar y0(const dimensionedScalar &ds)
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
Tensor< scalar > tensor
Definition symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
BarycentricTensor< scalar > barycentricTensor
A scalar version of the templated BarycentricTensor.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
#define registerInfoSwitch(Name, Type, SwitchVar)
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
3D tensor transformation operations.