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