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