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-2018 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::scalar Foam::particle::negativeSpaceDisplacementFactor = 1.01;
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  const bool isOwner = mesh_.faceOwner()[tetFacei_] == celli_;
156 
157  const label firstTetPtI = 1;
158  const label lastTetPtI = mesh_.faces()[tetFacei_].size() - 2;
159 
160  if (tetTriI == 1)
161  {
162  changeFace(tetTriI);
163  }
164  else if (tetTriI == 2)
165  {
166  if (isOwner)
167  {
168  if (tetPti_ == lastTetPtI)
169  {
170  changeFace(tetTriI);
171  }
172  else
173  {
174  reflect();
175  tetPti_ += 1;
176  }
177  }
178  else
179  {
180  if (tetPti_ == firstTetPtI)
181  {
182  changeFace(tetTriI);
183  }
184  else
185  {
186  reflect();
187  tetPti_ -= 1;
188  }
189  }
190  }
191  else if (tetTriI == 3)
192  {
193  if (isOwner)
194  {
195  if (tetPti_ == firstTetPtI)
196  {
197  changeFace(tetTriI);
198  }
199  else
200  {
201  reflect();
202  tetPti_ -= 1;
203  }
204  }
205  else
206  {
207  if (tetPti_ == lastTetPtI)
208  {
209  changeFace(tetTriI);
210  }
211  else
212  {
213  reflect();
214  tetPti_ += 1;
215  }
216  }
217  }
218  else
219  {
221  << "Changing tet without changing cell should only happen when the "
222  << "track is on triangle 1, 2 or 3."
223  << exit(FatalError);
224  }
225 }
226 
227 
228 void Foam::particle::changeFace(const label tetTriI)
229 {
230  // Get the old topology
231  const triFace triOldIs(currentTetIndices().faceTriIs(mesh_));
232 
233  // Get the shared edge and the pre-rotation
234  edge sharedEdge;
235  if (tetTriI == 1)
236  {
237  sharedEdge = edge(triOldIs[1], triOldIs[2]);
238  }
239  else if (tetTriI == 2)
240  {
241  sharedEdge = edge(triOldIs[2], triOldIs[0]);
242  }
243  else if (tetTriI == 3)
244  {
245  sharedEdge = edge(triOldIs[0], triOldIs[1]);
246  }
247  else
248  {
250  << "Changing face without changing cell should only happen when the"
251  << " track is on triangle 1, 2 or 3."
252  << exit(FatalError);
253 
254  sharedEdge = edge(-1, -1);
255  }
256 
257  // Find the face in the same cell that shares the edge, and the
258  // corresponding tetrahedra point
259  tetPti_ = -1;
260  forAll(mesh_.cells()[celli_], cellFaceI)
261  {
262  const label newFaceI = mesh_.cells()[celli_][cellFaceI];
263  const class face& newFace = mesh_.faces()[newFaceI];
264  const label newOwner = mesh_.faceOwner()[newFaceI];
265 
266  // Exclude the current face
267  if (tetFacei_ == newFaceI)
268  {
269  continue;
270  }
271 
272  // Loop over the edges, looking for the shared one. Note that we have to
273  // match the direction of the edge as well as the end points in order to
274  // avoid false positives when dealing with coincident ACMI faces.
275  const label edgeComp = newOwner == celli_ ? -1 : +1;
276  label edgeI = 0;
277  for
278  (
279  ;
280  edgeI < newFace.size()
281  && edge::compare(sharedEdge, newFace.faceEdge(edgeI)) != edgeComp;
282  ++ edgeI
283  );
284 
285  // If the face does not contain the edge, then move on to the next face
286  if (edgeI >= newFace.size())
287  {
288  continue;
289  }
290 
291  // Make the edge index relative to the base point
292  const label newBaseI = max(0, mesh_.tetBasePtIs()[newFaceI]);
293  edgeI = (edgeI - newBaseI + newFace.size()) % newFace.size();
294 
295  // If the edge is next the base point (i.e., the index is 0 or n - 1),
296  // then we swap it for the adjacent edge. This new edge is opposite the
297  // base point, and defines the tet with the original edge in it.
298  edgeI = min(max(1, edgeI), newFace.size() - 2);
299 
300  // Set the new face and tet point
301  tetFacei_ = newFaceI;
302  tetPti_ = edgeI;
303 
304  // Exit the loop now that the tet point has been found
305  break;
306  }
307 
308  if (tetPti_ == -1)
309  {
311  << "The search for an edge-connected face and tet-point failed."
312  << exit(FatalError);
313  }
314 
315  // Pre-rotation puts the shared edge opposite the base of the tetrahedron
316  if (sharedEdge.otherVertex(triOldIs[1]) == -1)
317  {
318  rotate(false);
319  }
320  else if (sharedEdge.otherVertex(triOldIs[2]) == -1)
321  {
322  rotate(true);
323  }
324 
325  // Get the new topology
326  const triFace triNewIs = currentTetIndices().faceTriIs(mesh_);
327 
328  // Reflect to account for the change of triangle orientation on the new face
329  reflect();
330 
331  // Post rotation puts the shared edge back in the correct location
332  if (sharedEdge.otherVertex(triNewIs[1]) == -1)
333  {
334  rotate(true);
335  }
336  else if (sharedEdge.otherVertex(triNewIs[2]) == -1)
337  {
338  rotate(false);
339  }
340 }
341 
342 
343 void Foam::particle::changeCell()
344 {
345  // Set the cell to be the one on the other side of the face
346  const label ownerCellI = mesh_.faceOwner()[tetFacei_];
347  const bool isOwner = celli_ == ownerCellI;
348  celli_ = isOwner ? mesh_.faceNeighbour()[tetFacei_] : ownerCellI;
349 
350  // Reflect to account for the change of triangle orientation in the new cell
351  reflect();
352 }
353 
354 
355 void Foam::particle::changeToMasterPatch()
356 {
357  label thisPatch = patch();
358 
359  forAll(mesh_.cells()[celli_], cellFaceI)
360  {
361  // Skip the current face and any internal faces
362  const label otherFaceI = mesh_.cells()[celli_][cellFaceI];
363  if (facei_ == otherFaceI || mesh_.isInternalFace(otherFaceI))
364  {
365  continue;
366  }
367 
368  // Compare the two faces. If they are the same, chose the one with the
369  // lower patch index. In the case of an ACMI-wall pair, this will be
370  // the ACMI, which is what we want.
371  const class face& thisFace = mesh_.faces()[facei_];
372  const class face& otherFace = mesh_.faces()[otherFaceI];
373  if (face::compare(thisFace, otherFace) != 0)
374  {
375  const label otherPatch =
376  mesh_.boundaryMesh().whichPatch(otherFaceI);
377  if (thisPatch > otherPatch)
378  {
379  facei_ = otherFaceI;
380  thisPatch = otherPatch;
381  }
382  }
383  }
384 
385  tetFacei_ = facei_;
386 }
387 
388 
389 void Foam::particle::locate
390 (
391  const vector& position,
392  const vector* direction,
393  label celli,
394  const bool boundaryFail,
395  const string boundaryMsg
396 )
397 {
398  // Find the cell, if it has not been given
399  if (celli < 0)
400  {
401  celli = mesh_.cellTree().findInside(position);
402  }
403  if (celli < 0)
404  {
406  << "Cell not found for particle position " << position << "."
407  << exit(FatalError);
408  }
409  celli_ = celli;
410 
411  // Track from the centre of the cell to the desired position
412  const vector displacement = position - mesh_.cellCentres()[celli_];
413 
414  // Loop all cell tets to find the one containing the position. Track
415  // through each tet from the cell centre. If a tet contains the position
416  // then the track will end with a single trackToTri.
417  const class cell& c = mesh_.cells()[celli_];
418  label minF = 1, minTetFacei = -1, minTetPti = -1;
419  forAll(c, cellTetFacei)
420  {
421  const class face& f = mesh_.faces()[c[cellTetFacei]];
422  for (label tetPti = 1; tetPti < f.size() - 1; ++ tetPti)
423  {
424  coordinates_ = barycentric(1, 0, 0, 0);
425  tetFacei_ = c[cellTetFacei];
426  tetPti_ = tetPti;
427  facei_ = -1;
428 
429  label tetTriI = -1;
430  const scalar f = trackToTri(displacement, 0, tetTriI);
431 
432  if (tetTriI == -1)
433  {
434  return;
435  }
436 
437  if (f < minF)
438  {
439  minF = f;
440  minTetFacei = tetFacei_;
441  minTetPti = tetPti_;
442  }
443  }
444  }
445 
446  // The particle must be (hopefully only slightly) outside the cell. Track
447  // into the tet which got the furthest.
448  coordinates_ = barycentric(1, 0, 0, 0);
449  tetFacei_ = minTetFacei;
450  tetPti_ = minTetPti;
451  facei_ = -1;
452 
453  track(displacement, 0);
454  if (!onFace())
455  {
456  return;
457  }
458 
459  // If we are here then we hit a boundary
460  if (boundaryFail)
461  {
462  FatalErrorInFunction << boundaryMsg << exit(FatalError);
463  }
464  else
465  {
466  // Re-do the track, but this time do the bit tangential to the
467  // direction/patch first. This gets us as close as possible to the
468  // original path/position.
469 
470  if (direction == nullptr)
471  {
472  const polyPatch& p = mesh_.boundaryMesh()[patch()];
473  direction = &p.faceNormals()[p.whichFace(facei_)];
474  }
475 
476  const vector n = *direction/mag(*direction);
477  const vector sN = (displacement & n)*n;
478  const vector sT = displacement - sN;
479 
480  coordinates_ = barycentric(1, 0, 0, 0);
481  celli_ = celli;
482  tetFacei_ = minTetFacei;
483  tetPti_ = minTetPti;
484  facei_ = -1;
485 
486  track(sT, 0);
487  track(sN, 0);
488 
489  static label nWarnings = 0;
490  static const label maxNWarnings = 100;
491  if (nWarnings < maxNWarnings)
492  {
493  WarningInFunction << boundaryMsg << endl;
494  ++ nWarnings;
495  }
496  if (nWarnings == maxNWarnings)
497  {
499  << "Suppressing any further warnings about particles being "
500  << "located outside of the mesh." << endl;
501  ++ nWarnings;
502  }
503  }
504 }
505 
506 
507 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
508 
510 (
511  const polyMesh& mesh,
512  const barycentric& coordinates,
513  const label celli,
514  const label tetFacei,
515  const label tetPti
516 )
517 :
518  mesh_(mesh),
519  coordinates_(coordinates),
520  celli_(celli),
521  tetFacei_(tetFacei),
522  tetPti_(tetPti),
523  facei_(-1),
524  stepFraction_(0.0),
525  origProc_(Pstream::myProcNo()),
526  origId_(getNewParticleID())
527 {}
528 
529 
531 (
532  const polyMesh& mesh,
533  const vector& position,
534  const label celli
535 )
536 :
537  mesh_(mesh),
538  coordinates_(- vGreat, - vGreat, - vGreat, - vGreat),
539  celli_(celli),
540  tetFacei_(-1),
541  tetPti_(-1),
542  facei_(-1),
543  stepFraction_(0.0),
544  origProc_(Pstream::myProcNo()),
545  origId_(getNewParticleID())
546 {
547  locate
548  (
549  position,
550  nullptr,
551  celli,
552  false,
553  "Particle initialised with a location outside of the mesh."
554  );
555 }
556 
557 
559 :
560  mesh_(p.mesh_),
561  coordinates_(p.coordinates_),
562  celli_(p.celli_),
563  tetFacei_(p.tetFacei_),
564  tetPti_(p.tetPti_),
565  facei_(p.facei_),
566  stepFraction_(p.stepFraction_),
567  origProc_(p.origProc_),
568  origId_(p.origId_)
569 {}
570 
571 
573 :
574  mesh_(mesh),
575  coordinates_(p.coordinates_),
576  celli_(p.celli_),
577  tetFacei_(p.tetFacei_),
578  tetPti_(p.tetPti_),
579  facei_(p.facei_),
580  stepFraction_(p.stepFraction_),
581  origProc_(p.origProc_),
582  origId_(p.origId_)
583 {}
584 
585 
586 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
587 
588 Foam::scalar Foam::particle::track
589 (
590  const vector& displacement,
591  const scalar fraction
592 )
593 {
594  scalar f = trackToFace(displacement, fraction);
595 
596  while (onInternalFace())
597  {
598  changeCell();
599 
600  f *= trackToFace(f*displacement, f*fraction);
601  }
602 
603  return f;
604 }
605 
606 
607 Foam::scalar Foam::particle::trackToCell
608 (
609  const vector& displacement,
610  const scalar fraction
611 )
612 {
613  const scalar f = trackToFace(displacement, fraction);
614 
615  if (onInternalFace())
616  {
617  changeCell();
618  }
619 
620  return f;
621 }
622 
623 
624 Foam::scalar Foam::particle::trackToFace
625 (
626  const vector& displacement,
627  const scalar fraction
628 )
629 {
630  scalar f = 1;
631 
632  label tetTriI = onFace() ? 0 : -1;
633 
634  facei_ = -1;
635 
636  while (true)
637  {
638  f *= trackToTri(f*displacement, f*fraction, tetTriI);
639 
640  if (tetTriI == -1)
641  {
642  // The track has completed within the current tet
643  return 0;
644  }
645  else if (tetTriI == 0)
646  {
647  // The track has hit a face, so set the current face and return
648  facei_ = tetFacei_;
649  return f;
650  }
651  else
652  {
653  // Move to the next tet and continue the track
654  changeTet(tetTriI);
655  }
656  }
657 }
658 
659 
661 (
662  const vector& displacement,
663  const scalar fraction,
664  label& tetTriI
665 )
666 {
667  const vector x0 = position();
668  const vector x1 = displacement;
669  const barycentric y0 = coordinates_;
670 
671  if (debug)
672  {
673  Info<< "Particle " << origId() << endl << "Tracking from " << x0
674  << " along " << x1 << " to " << x0 + x1 << endl;
675  }
676 
677  // Get the tet geometry
678  vector centre;
679  scalar detA;
681  stationaryTetReverseTransform(centre, detA, T);
682 
683  if (debug)
684  {
685  vector o, b, v1, v2;
686  stationaryTetGeometry(o, b, v1, v2);
687  Info<< "Tet points o=" << o << ", b=" << b
688  << ", v1=" << v1 << ", v2=" << v2 << endl
689  << "Tet determinant = " << detA << endl
690  << "Start local coordinates = " << y0 << endl;
691  }
692 
693  // Get the factor by which the displacement is increased
694  const scalar f = detA >= 0 ? 1 : negativeSpaceDisplacementFactor;
695 
696  // Calculate the local tracking displacement
697  barycentric Tx1(f*x1 & T);
698 
699  if (debug)
700  {
701  Info<< "Local displacement = " << Tx1 << "/" << detA << endl;
702  }
703 
704  // Calculate the hit fraction
705  label iH = -1;
706  scalar muH = std::isnormal(detA) && detA <= 0 ? vGreat : 1/detA;
707  for (label i = 0; i < 4; ++ i)
708  {
709  if (std::isnormal(Tx1[i]) && Tx1[i] < 0)
710  {
711  scalar mu = - y0[i]/Tx1[i];
712 
713  if (debug)
714  {
715  Info<< "Hit on tet face " << i << " at local coordinate "
716  << y0 + mu*Tx1 << ", " << mu*detA*100 << "% of the "
717  << "way along the track" << endl;
718  }
719 
720  if (0 <= mu && mu < muH)
721  {
722  iH = i;
723  muH = mu;
724  }
725  }
726  }
727 
728  // Set the new coordinates
729  barycentric yH = y0 + muH*Tx1;
730 
731  // Clamp to zero any negative coordinates generated by round-off error
732  for (label i = 0; i < 4; ++ i)
733  {
734  yH.replace(i, i == iH ? 0 : max(0, yH[i]));
735  }
736 
737  // Re-normalise if within the tet
738  if (iH == -1)
739  {
740  yH /= cmptSum(yH);
741  }
742 
743  // Set the new position and hit index
744  coordinates_ = yH;
745  tetTriI = iH;
746 
747  if (debug)
748  {
749  if (iH != -1)
750  {
751  Info<< "Track hit tet face " << iH << " first" << endl;
752  }
753  else
754  {
755  Info<< "Track hit no tet faces" << endl;
756  }
757  Info<< "End local coordinates = " << yH << endl
758  << "End global coordinates = " << position() << endl
759  << "Tracking displacement = " << position() - x0 << endl
760  << muH*detA*100 << "% of the step from " << stepFraction_ << " to "
761  << stepFraction_ + fraction << " completed" << endl << endl;
762  }
763 
764  // Set the proportion of the track that has been completed
765  stepFraction_ += fraction*muH*detA;
766 
767  return iH != -1 ? 1 - muH*detA : 0;
768 }
769 
770 
772 (
773  const vector& displacement,
774  const scalar fraction,
775  label& tetTriI
776 )
777 {
778  const vector x0 = position();
779  const vector x1 = displacement;
780  const barycentric y0 = coordinates_;
781 
782  if (debug)
783  {
784  Info<< "Particle " << origId() << endl << "Tracking from " << x0
785  << " along " << x1 << " to " << x0 + x1 << endl;
786  }
787 
788  // Get the tet geometry
789  Pair<vector> centre;
792  movingTetReverseTransform(fraction, centre, detA, T);
793 
794  if (debug)
795  {
796  Pair<vector> o, b, v1, v2;
797  movingTetGeometry(fraction, o, b, v1, v2);
798  Info<< "Tet points o=" << o[0] << ", b=" << b[0]
799  << ", v1=" << v1[0] << ", v2=" << v2[0] << endl
800  << "Tet determinant = " << detA[0] << endl
801  << "Start local coordinates = " << y0[0] << endl;
802  }
803 
804  // Get the factor by which the displacement is increased
805  const scalar f = detA[0] >= 0 ? 1 : negativeSpaceDisplacementFactor;
806 
807  // Get the relative global position
808  const vector x0Rel = x0 - centre[0];
809  const vector x1Rel = f*x1 - centre[1];
810 
811  // Form the determinant and hit equations
812  cubicEqn detAEqn(sqr(detA[0])*detA[3], detA[0]*detA[2], detA[1], 1);
813  const barycentric yC(1, 0, 0, 0);
814  const barycentric hitEqnA =
815  ((x1Rel & T[2]) + detA[3]*yC)*sqr(detA[0]);
816  const barycentric hitEqnB =
817  ((x1Rel & T[1]) + (x0Rel & T[2]) + detA[2]*yC)*detA[0];
818  const barycentric hitEqnC =
819  ((x1Rel & T[0]) + (x0Rel & T[1]) + detA[1]*yC);
820  const barycentric hitEqnD = y0;
821  FixedList<cubicEqn, 4> hitEqn;
822  forAll(hitEqn, i)
823  {
824  hitEqn[i] = cubicEqn(hitEqnA[i], hitEqnB[i], hitEqnC[i], hitEqnD[i]);
825  }
826 
827  if (debug)
828  {
829  for (label i = 0; i < 4; ++ i)
830  {
831  Info<< (i ? " " : "Hit equation ") << i << " = "
832  << hitEqn[i] << endl;
833  }
834  Info<< " DetA equation = " << detA << endl;
835  }
836 
837  // Calculate the hit fraction
838  label iH = -1;
839  scalar muH = std::isnormal(detA[0]) && detA[0] <= 0 ? vGreat : 1/detA[0];
840  for (label i = 0; i < 4; ++ i)
841  {
842  const Roots<3> mu = hitEqn[i].roots();
843 
844  for (label j = 0; j < 3; ++ j)
845  {
846  if (mu.type(j) == roots::real && hitEqn[i].derivative(mu[j]) < 0)
847  {
848  if (debug)
849  {
850  const barycentric yH
851  (
852  hitEqn[0].value(mu[j]),
853  hitEqn[1].value(mu[j]),
854  hitEqn[2].value(mu[j]),
855  hitEqn[3].value(mu[j])
856  );
857  const scalar detAH = detAEqn.value(mu[j]);
858 
859  Info<< "Hit on tet face " << i << " at local coordinate "
860  << yH/detAH << ", " << mu[j]*detA[0]*100 << "% of the "
861  << "way along the track" << endl;
862  }
863 
864  if (0 <= mu[j] && mu[j] < muH)
865  {
866  iH = i;
867  muH = mu[j];
868  }
869  }
870  }
871  }
872 
873  // Set the new coordinates
874  barycentric yH
875  (
876  hitEqn[0].value(muH),
877  hitEqn[1].value(muH),
878  hitEqn[2].value(muH),
879  hitEqn[3].value(muH)
880  );
881  // !!! <-- This fails if the tet collapses onto the particle, as detA tends
882  // to zero at the hit. In this instance, we can differentiate the hit and
883  // detA polynomials to find a limiting location, but this will not be on a
884  // triangle. We will then need to do a second track through the degenerate
885  // tet to find the final hit position. This second track is over zero
886  // distance and therefore can be of the static mesh type. This has not yet
887  // been implemented.
888  const scalar detAH = detAEqn.value(muH);
889  if (!std::isnormal(detAH))
890  {
892  << "A moving tet collapsed onto a particle. This is not supported. "
893  << "The mesh is too poor, or the motion too severe, for particle "
894  << "tracking to function." << exit(FatalError);
895  }
896  yH /= detAH;
897 
898  // Clamp to zero any negative coordinates generated by round-off error
899  for (label i = 0; i < 4; ++ i)
900  {
901  yH.replace(i, i == iH ? 0 : max(0, yH[i]));
902  }
903 
904  // Re-normalise if within the tet
905  if (iH == -1)
906  {
907  yH /= cmptSum(yH);
908  }
909 
910  // Set the new position and hit index
911  coordinates_ = yH;
912  tetTriI = iH;
913 
914  // Set the proportion of the track that has been completed
915  stepFraction_ += fraction*muH*detA[0];
916 
917  if (debug)
918  {
919  if (iH != -1)
920  {
921  Info<< "Track hit tet face " << iH << " first" << endl;
922  }
923  else
924  {
925  Info<< "Track hit no tet faces" << endl;
926  }
927  Info<< "End local coordinates = " << yH << endl
928  << "End global coordinates = " << position() << endl
929  << "Tracking displacement = " << position() - x0 << endl
930  << muH*detA[0]*100 << "% of the step from " << stepFraction_
931  << " to " << stepFraction_ + fraction << " completed" << endl
932  << endl;
933  }
934 
935  return iH != -1 ? 1 - muH*detA[0] : 0;
936 }
937 
938 
939 Foam::scalar Foam::particle::trackToTri
940 (
941  const vector& displacement,
942  const scalar fraction,
943  label& tetTriI
944 )
945 {
946  if (mesh_.moving())
947  {
948  return trackToMovingTri(displacement, fraction, tetTriI);
949  }
950  else
951  {
952  return trackToStationaryTri(displacement, fraction, tetTriI);
953  }
954 }
955 
956 
958 {
959  if (cmptMin(mesh_.geometricD()) == -1)
960  {
961  vector pos = position(), posC = pos;
963  return pos - posC;
964  }
965  else
966  {
967  return vector::zero;
968  }
969 }
970 
971 
973 {}
974 
975 
977 {}
978 
979 
981 {
982  // Convert the face index to be local to the processor patch
983  facei_ = mesh_.boundaryMesh()[patch()].whichFace(facei_);
984 }
985 
986 
988 (
989  const label patchi,
990  trackingData& td
991 )
992 {
993  const coupledPolyPatch& ppp =
994  refCast<const coupledPolyPatch>(mesh_.boundaryMesh()[patchi]);
995 
996  if (!ppp.parallel())
997  {
998  const tensor& T =
999  (
1000  ppp.forwardT().size() == 1
1001  ? ppp.forwardT()[0]
1002  : ppp.forwardT()[facei_]
1003  );
1005  }
1006  else if (ppp.separated())
1007  {
1008  const vector& s =
1009  (
1010  (ppp.separation().size() == 1)
1011  ? ppp.separation()[0]
1012  : ppp.separation()[facei_]
1013  );
1014  transformProperties(-s);
1015  }
1016 
1017  // Set the topology
1018  celli_ = ppp.faceCells()[facei_];
1019  facei_ += ppp.start();
1020  tetFacei_ = facei_;
1021  // Faces either side of a coupled patch are numbered in opposite directions
1022  // as their normals both point away from their connected cells. The tet
1023  // point therefore counts in the opposite direction from the base point.
1024  tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;
1025 
1026  // Reflect to account for the change of triangle orientation in the new cell
1027  reflect();
1028 
1029  // Note that the position does not need transforming explicitly. The face-
1030  // triangle on the receive patch is the transformation of the one on the
1031  // send patch, so whilst the barycentric coordinates remain the same, the
1032  // change of triangle implicitly transforms the position.
1033 }
1034 
1035 
1038  const vectorTensorTransform& transform
1039 )
1040 {
1041  // Get the transformed position
1042  const vector pos = transform.invTransformPosition(position());
1043 
1044  // Break the topology
1045  celli_ = -1;
1046  tetFacei_ = -1;
1047  tetPti_ = -1;
1048  facei_ = -1;
1049 
1050  // Store the position in the barycentric data
1051  coordinates_ = barycentric(1 - cmptSum(pos), pos.x(), pos.y(), pos.z());
1052 
1053  // Transform the properties
1054  transformProperties(- transform.t());
1055  if (transform.hasR())
1056  {
1057  transformProperties(transform.R().T());
1058  }
1059 }
1060 
1061 
1063 {
1064  // Get the position from the barycentric data
1065  const vector pos(coordinates_.b(), coordinates_.c(), coordinates_.d());
1066 
1067  // Create some arbitrary topology for the supplied cell
1068  celli_ = celli;
1069  tetFacei_ = mesh_.cells()[celli_][0];
1070  tetPti_ = 1;
1071  facei_ = -1;
1072 
1073  // Get the reverse transform and directly set the coordinates from the
1074  // position. This isn't likely to be correct; the particle is probably not
1075  // in this tet. It will, however, generate the correct vector when the
1076  // position method is called. A referred particle should never be tracked,
1077  // so this approximate topology is good enough. By using the nearby cell we
1078  // minimize the error associated with the incorrect topology.
1079  coordinates_ = barycentric(1, 0, 0, 0);
1080  if (mesh_.moving())
1081  {
1082  Pair<vector> centre;
1083  FixedList<scalar, 4> detA;
1085  movingTetReverseTransform(0, centre, detA, T);
1086  coordinates_ += (pos - centre[0]) & T[0]/detA[0];
1087  }
1088  else
1089  {
1090  vector centre;
1091  scalar detA;
1093  stationaryTetReverseTransform(centre, detA, T);
1094  coordinates_ += (pos - centre) & T/detA;
1095  }
1096 }
1097 
1098 
1101  const polyMesh& procMesh,
1102  const label procCell,
1103  const label procTetFace
1104 ) const
1105 {
1106  // The tet point on the procMesh differs from the current tet point if the
1107  // mesh and procMesh faces are of differing orientation. The change is the
1108  // same as in particle::correctAfterParallelTransfer.
1109 
1110  if
1111  (
1112  (mesh_.faceOwner()[tetFacei_] == celli_)
1113  == (procMesh.faceOwner()[procTetFace] == procCell)
1114  )
1115  {
1116  return tetPti_;
1117  }
1118  else
1119  {
1120  return procMesh.faces()[procTetFace].size() - 1 - tetPti_;
1121  }
1122 }
1123 
1124 
1127  const vector& position,
1128  const mapPolyMesh& mapper
1129 )
1130 {
1131  locate
1132  (
1133  position,
1134  nullptr,
1135  mapper.reverseCellMap()[celli_],
1136  true,
1137  "Particle mapped to a location outside of the mesh."
1138  );
1139 }
1140 
1141 
1142 // * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
1143 
1144 bool Foam::operator==(const particle& pA, const particle& pB)
1145 {
1146  return (pA.origProc() == pB.origProc() && pA.origId() == pB.origId());
1147 }
1148 
1149 
1150 bool Foam::operator!=(const particle& pA, const particle& pB)
1151 {
1152  return !(pA == pB);
1153 }
1154 
1155 
1156 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
label origProc() const
Return the originating processor ID.
Definition: particleI.H:178
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition: particle.C:988
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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 deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:957
virtual bool separated() const
Are the planes separated.
bool moving() const
Is mesh moving.
Definition: polyMesh.H:502
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:1100
Vector-tensor class used to perform translations and rotations in 3D space.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const Cmpt & b() const
Definition: BarycentricI.H:66
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:625
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Tensor< Cmpt > T() const
Return transpose.
Definition: TensorI.H:321
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:940
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:427
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:500
vector invTransformPosition(const vector &v) const
Inverse transform the given position.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
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:68
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 void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:972
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:608
virtual const vectorField & separation() const
If the planes are separated the separation vector.
void prepareForInteractionListReferral(const vectorTensorTransform &transform)
Break the topology and store the particle position so that the.
Definition: particle.C:1037
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual const tensorField & forwardT() const
Return face transformation tensor.
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:807
virtual bool parallel() const
Are the cyclic planes parallel.
dimensionedScalar pos(const dimensionedScalar &ds)
Cubic equation of the form a*x^3 + b*x^2 + c*x + d = 0.
Definition: cubicEqn.H:49
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
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:315
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:510
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:1062
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1041
const Cmpt & d() const
Definition: BarycentricI.H:80
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
Definition: particle.C:980
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1028
const Cmpt & x() const
Definition: VectorI.H:75
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:190
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1126
const Cmpt & c() const
Definition: BarycentricI.H:73
bool onFace() const
Is the particle on a face?
Definition: particleI.H:259
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:526
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.
const dimensionedScalar mu
Atomic mass unit.
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:772
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
void type(const direction i, const roots::type t)
Set the type of the i-th root.
Definition: RootsI.H:120
const dimensionedScalar c
Speed of light in a vacuum.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
static label particleCount_
Cumulative particle counter - used to provode unique ID.
Definition: particle.H:349
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:265
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:661
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:589
bool operator!=(const particle &, const particle &)
Definition: particle.C:1150
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:283
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:277
Namespace for OpenFOAM.