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