particle.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2017 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::stationaryTetGeometry
46 (
47  vector& centre,
48  vector& base,
49  vector& vertex1,
50  vector& vertex2
51 ) const
52 {
53  const triFace triIs(currentTetIndices().faceTriIs(mesh_));
54  const vectorField& ccs = mesh_.cellCentres();
55  const pointField& pts = mesh_.points();
56 
57  centre = ccs[celli_];
58  base = pts[triIs[0]];
59  vertex1 = pts[triIs[1]];
60  vertex2 = pts[triIs[2]];
61 }
62 
63 
64 Foam::barycentricTensor Foam::particle::stationaryTetTransform() const
65 {
66  vector centre, base, vertex1, vertex2;
67  stationaryTetGeometry(centre, base, vertex1, vertex2);
68 
69  return barycentricTensor(centre, base, vertex1, vertex2);
70 }
71 
72 
73 void Foam::particle::stationaryTetReverseTransform
74 (
75  vector& centre,
76  scalar& detA,
78 ) const
79 {
80  barycentricTensor A = stationaryTetTransform();
81 
82  const vector ab = A.b() - A.a();
83  const vector ac = A.c() - A.a();
84  const vector ad = A.d() - A.a();
85  const vector bc = A.c() - A.b();
86  const vector bd = A.d() - A.b();
87 
88  centre = A.a();
89 
90  detA = ab & (ac ^ ad);
91 
93  (
94  bd ^ bc,
95  ac ^ ad,
96  ad ^ ab,
97  ab ^ ac
98  );
99 }
100 
101 
102 void Foam::particle::movingTetGeometry
103 (
104  const scalar fraction,
105  Pair<vector>& centre,
106  Pair<vector>& base,
107  Pair<vector>& vertex1,
108  Pair<vector>& vertex2
109 ) const
110 {
111  const triFace triIs(currentTetIndices().faceTriIs(mesh_));
112  const pointField& ptsOld = mesh_.oldPoints();
113  const pointField& ptsNew = mesh_.points();
114 
115  // !!! <-- We would be better off using mesh_.cellCentres() here. However,
116  // we need to put a mesh_.oldCellCentres() method in for this to work. The
117  // values obtained from the mesh and those obtained from the cell do not
118  // necessarily match. See mantis #1993.
119  const vector ccOld = mesh_.cells()[celli_].centre(ptsOld, mesh_.faces());
120  const vector ccNew = mesh_.cells()[celli_].centre(ptsNew, mesh_.faces());
121 
122  // Old and new points and cell centres are not sub-cycled. If we are sub-
123  // cycling, then we have to account for the timestep change here by
124  // modifying the fractions that we take of the old and new geometry.
125  const Pair<scalar> s = stepFractionSpan();
126  const scalar f0 = s[0] + stepFraction_*s[1], f1 = fraction*s[1];
127 
128  centre[0] = ccOld + f0*(ccNew - ccOld);
129  base[0] = ptsOld[triIs[0]] + f0*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
130  vertex1[0] = ptsOld[triIs[1]] + f0*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
131  vertex2[0] = ptsOld[triIs[2]] + f0*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
132 
133  centre[1] = f1*(ccNew - ccOld);
134  base[1] = f1*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
135  vertex1[1] = f1*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
136  vertex2[1] = f1*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
137 }
138 
139 
140 Foam::Pair<Foam::barycentricTensor> Foam::particle::movingTetTransform
141 (
142  const scalar fraction
143 ) const
144 {
145  Pair<vector> centre, base, vertex1, vertex2;
146  movingTetGeometry(fraction, centre, base, vertex1, vertex2);
147 
148  return
149  Pair<barycentricTensor>
150  (
151  barycentricTensor(centre[0], base[0], vertex1[0], vertex2[0]),
152  barycentricTensor(centre[1], base[1], vertex1[1], vertex2[1])
153  );
154 }
155 
156 
157 void Foam::particle::movingTetReverseTransform
158 (
159  const scalar fraction,
160  Pair<vector>& centre,
161  FixedList<scalar, 4>& detA,
162  FixedList<barycentricTensor, 3>& T
163 ) const
164 {
165  Pair<barycentricTensor> A = movingTetTransform(fraction);
166 
167  const Pair<vector> ab(A[0].b() - A[0].a(), A[1].b() - A[1].a());
168  const Pair<vector> ac(A[0].c() - A[0].a(), A[1].c() - A[1].a());
169  const Pair<vector> ad(A[0].d() - A[0].a(), A[1].d() - A[1].a());
170  const Pair<vector> bc(A[0].c() - A[0].b(), A[1].c() - A[1].b());
171  const Pair<vector> bd(A[0].d() - A[0].b(), A[1].d() - A[1].b());
172 
173  centre[0] = A[0].a();
174  centre[1] = A[1].a();
175 
176  detA[0] = ab[0] & (ac[0] ^ ad[0]);
177  detA[1] =
178  (ab[1] & (ac[0] ^ ad[0]))
179  + (ab[0] & (ac[1] ^ ad[0]))
180  + (ab[0] & (ac[0] ^ ad[1]));
181  detA[2] =
182  (ab[0] & (ac[1] ^ ad[1]))
183  + (ab[1] & (ac[0] ^ ad[1]))
184  + (ab[1] & (ac[1] ^ ad[0]));
185  detA[3] = ab[1] & (ac[1] ^ ad[1]);
186 
187  T[0] = barycentricTensor
188  (
189  bd[0] ^ bc[0],
190  ac[0] ^ ad[0],
191  ad[0] ^ ab[0],
192  ab[0] ^ ac[0]
193  );
194  T[1] = barycentricTensor
195  (
196  (bd[0] ^ bc[1]) + (bd[1] ^ bc[0]),
197  (ac[0] ^ ad[1]) + (ac[1] ^ ad[0]),
198  (ad[0] ^ ab[1]) + (ad[1] ^ ab[0]),
199  (ab[0] ^ ac[1]) + (ab[1] ^ ac[0])
200  );
201  T[2] = barycentricTensor
202  (
203  bd[1] ^ bc[1],
204  ac[1] ^ ad[1],
205  ad[1] ^ ab[1],
206  ab[1] ^ ac[1]
207  );
208 }
209 
210 
211 void Foam::particle::reflect()
212 {
213  Swap(coordinates_.c(), coordinates_.d());
214 }
215 
216 
217 void Foam::particle::rotate(const bool reverse)
218 {
219  if (!reverse)
220  {
221  scalar temp = coordinates_.b();
222  coordinates_.b() = coordinates_.c();
223  coordinates_.c() = coordinates_.d();
224  coordinates_.d() = temp;
225  }
226  else
227  {
228  scalar temp = coordinates_.d();
229  coordinates_.d() = coordinates_.c();
230  coordinates_.c() = coordinates_.b();
231  coordinates_.b() = temp;
232  }
233 }
234 
235 
236 void Foam::particle::changeTet(const label tetTriI)
237 {
238  const bool isOwner = mesh_.faceOwner()[tetFacei_] == celli_;
239 
240  const label firstTetPtI = 1;
241  const label lastTetPtI = mesh_.faces()[tetFacei_].size() - 2;
242 
243  if (tetTriI == 1)
244  {
245  changeFace(tetTriI);
246  }
247  else if (tetTriI == 2)
248  {
249  if (isOwner)
250  {
251  if (tetPti_ == lastTetPtI)
252  {
253  changeFace(tetTriI);
254  }
255  else
256  {
257  reflect();
258  tetPti_ += 1;
259  }
260  }
261  else
262  {
263  if (tetPti_ == firstTetPtI)
264  {
265  changeFace(tetTriI);
266  }
267  else
268  {
269  reflect();
270  tetPti_ -= 1;
271  }
272  }
273  }
274  else if (tetTriI == 3)
275  {
276  if (isOwner)
277  {
278  if (tetPti_ == firstTetPtI)
279  {
280  changeFace(tetTriI);
281  }
282  else
283  {
284  reflect();
285  tetPti_ -= 1;
286  }
287  }
288  else
289  {
290  if (tetPti_ == lastTetPtI)
291  {
292  changeFace(tetTriI);
293  }
294  else
295  {
296  reflect();
297  tetPti_ += 1;
298  }
299  }
300  }
301  else
302  {
304  << "Changing tet without changing cell should only happen when the "
305  << "track is on triangle 1, 2 or 3."
306  << exit(FatalError);
307  }
308 }
309 
310 
311 void Foam::particle::changeFace(const label tetTriI)
312 {
313  // Get the old topology
314  const triFace triOldIs(currentTetIndices().faceTriIs(mesh_));
315 
316  // Get the shared edge and the pre-rotation
317  edge sharedEdge;
318  if (tetTriI == 1)
319  {
320  sharedEdge = edge(triOldIs[1], triOldIs[2]);
321  }
322  else if (tetTriI == 2)
323  {
324  sharedEdge = edge(triOldIs[2], triOldIs[0]);
325  }
326  else if (tetTriI == 3)
327  {
328  sharedEdge = edge(triOldIs[0], triOldIs[1]);
329  }
330  else
331  {
333  << "Changing face without changing cell should only happen when the"
334  << " track is on triangle 1, 2 or 3."
335  << exit(FatalError);
336 
337  sharedEdge = edge(-1, -1);
338  }
339 
340  // Find the face in the same cell that shares the edge, and the
341  // corresponding tetrahedra point
342  tetPti_ = -1;
343  forAll(mesh_.cells()[celli_], cellFaceI)
344  {
345  const label newFaceI = mesh_.cells()[celli_][cellFaceI];
346  const class face& newFace = mesh_.faces()[newFaceI];
347 
348  // Exclude the current face
349  if (tetFacei_ == newFaceI)
350  {
351  continue;
352  }
353 
354  // Loop over the edges, looking for the shared one
355  label edgeComp = 0;
356  label edgeI = 0;
357  for (; edgeI < newFace.size(); ++ edgeI)
358  {
359  edgeComp = edge::compare(sharedEdge, newFace.faceEdge(edgeI));
360 
361  if (edgeComp != 0)
362  {
363  break;
364  }
365  }
366 
367  // If the face does not contain the edge, then move on to the next face
368  if (edgeComp == 0)
369  {
370  continue;
371  }
372 
373  // Make the edge index relative to the base point
374  const label newBaseI = max(0, mesh_.tetBasePtIs()[newFaceI]);
375  edgeI = (edgeI - newBaseI + newFace.size()) % newFace.size();
376 
377  // If the edge is next the base point (i.e., the index is 0 or n - 1),
378  // then we swap it for the adjacent edge. This new edge is opposite the
379  // base point, and defines the tet with the original edge in it.
380  edgeI = min(max(1, edgeI), newFace.size() - 2);
381 
382  // Set the new face and tet point
383  tetFacei_ = newFaceI;
384  tetPti_ = edgeI;
385 
386  // Exit the loop now that the the tet point has been found
387  break;
388  }
389 
390  if (tetPti_ == -1)
391  {
393  << "The search for an edge-connected face and tet-point failed."
394  << exit(FatalError);
395  }
396 
397  // Pre-rotation puts the shared edge opposite the base of the tetrahedron
398  if (sharedEdge.otherVertex(triOldIs[1]) == -1)
399  {
400  rotate(false);
401  }
402  else if (sharedEdge.otherVertex(triOldIs[2]) == -1)
403  {
404  rotate(true);
405  }
406 
407  // Get the new topology
408  const triFace triNewIs = currentTetIndices().faceTriIs(mesh_);
409 
410  // Reflect to account for the change of triangle orientation on the new face
411  reflect();
412 
413  // Post rotation puts the shared edge back in the correct location
414  if (sharedEdge.otherVertex(triNewIs[1]) == -1)
415  {
416  rotate(true);
417  }
418  else if (sharedEdge.otherVertex(triNewIs[2]) == -1)
419  {
420  rotate(false);
421  }
422 }
423 
424 
425 void Foam::particle::changeCell()
426 {
427  // Set the cell to be the one on the other side of the face
428  const label ownerCellI = mesh_.faceOwner()[tetFacei_];
429  const bool isOwner = celli_ == ownerCellI;
430  celli_ = isOwner ? mesh_.faceNeighbour()[tetFacei_] : ownerCellI;
431 
432  // Reflect to account for the change of triangle orientation in the new cell
433  reflect();
434 }
435 
436 
437 void Foam::particle::locate
438 (
439  const vector& position,
440  const vector* direction,
441  const label celli,
442  const bool boundaryFail,
443  const string boundaryMsg
444 )
445 {
446  celli_ = celli;
447 
448  // Find the cell, if it has not been given
449  if (celli_ < 0)
450  {
451  celli_ = mesh_.cellTree().findInside(position);
452  }
453  if (celli_ < 0)
454  {
456  << "Cell not found for particle position " << position << "."
457  << exit(FatalError);
458  }
459 
460  // Put the particle at the cell centre and in a random tet
461  coordinates_ = barycentric(1, 0, 0, 0);
462  tetFacei_ = mesh_.cells()[celli_][0];
463  tetPti_ = 1;
464  facei_ = -1;
465 
466  // Track to the injection point
467  track(position - mesh_.cellCentres()[celli_], 0);
468  if (!onFace())
469  {
470  return;
471  }
472 
473  // We hit a boundary ...
474  if (boundaryFail)
475  {
476  FatalErrorInFunction << boundaryMsg << exit(FatalError);
477  }
478  else
479  {
480  // Re-do the track, but this time do the bit tangential to the
481  // direction/patch first. This gets us as close as possible to the
482  // original path/position.
483 
484  if (direction == nullptr)
485  {
486  const polyPatch& p = mesh_.boundaryMesh()[patch()];
487  direction = &p.faceNormals()[p.whichFace(facei_)];
488  }
489 
490  const vector n = *direction/mag(*direction);
491  const vector s = position - mesh_.cellCentres()[celli_];
492  const vector sN = (s & n)*n;
493  const vector sT = s - sN;
494 
495  coordinates_ = barycentric(1, 0, 0, 0);
496  tetFacei_ = mesh_.cells()[celli_][0];
497  tetPti_ = 1;
498  facei_ = -1;
499 
500  track(sT, 0);
501  track(sN, 0);
502 
503  static label nWarnings = 0;
504  static const label maxNWarnings = 100;
505  if (nWarnings < maxNWarnings)
506  {
507  WarningInFunction << boundaryMsg << endl;
508  ++ nWarnings;
509  }
510  if (nWarnings == maxNWarnings)
511  {
513  << "Suppressing any further warnings about particles being "
514  << "located outside of the mesh." << endl;
515  ++ nWarnings;
516  }
517  }
518 }
519 
520 
521 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
522 
524 (
525  const polyMesh& mesh,
526  const barycentric& coordinates,
527  const label celli,
528  const label tetFacei,
529  const label tetPti
530 )
531 :
532  mesh_(mesh),
533  coordinates_(coordinates),
534  celli_(celli),
535  tetFacei_(tetFacei),
536  tetPti_(tetPti),
537  facei_(-1),
538  stepFraction_(0.0),
539  origProc_(Pstream::myProcNo()),
540  origId_(getNewParticleID())
541 {}
542 
543 
545 (
546  const polyMesh& mesh,
547  const vector& position,
548  const label celli
549 )
550 :
551  mesh_(mesh),
552  coordinates_(- VGREAT, - VGREAT, - VGREAT, - VGREAT),
553  celli_(celli),
554  tetFacei_(-1),
555  tetPti_(-1),
556  facei_(-1),
557  stepFraction_(0.0),
558  origProc_(Pstream::myProcNo()),
559  origId_(getNewParticleID())
560 {
561  locate
562  (
563  position,
564  nullptr,
565  celli,
566  false,
567  "Particle initialised with a location outside of the mesh."
568  );
569 }
570 
571 
573 :
574  mesh_(p.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 
587 :
588  mesh_(mesh),
589  coordinates_(p.coordinates_),
590  celli_(p.celli_),
591  tetFacei_(p.tetFacei_),
592  tetPti_(p.tetPti_),
593  facei_(p.facei_),
594  stepFraction_(p.stepFraction_),
595  origProc_(p.origProc_),
596  origId_(p.origId_)
597 {}
598 
599 
600 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
601 
602 Foam::scalar Foam::particle::track
603 (
604  const vector& displacement,
605  const scalar fraction
606 )
607 {
608  scalar f = trackToFace(displacement, fraction);
609 
610  while (onInternalFace())
611  {
612  changeCell();
613 
614  f *= trackToFace(f*displacement, f*fraction);
615  }
616 
617  return f;
618 }
619 
620 
621 Foam::scalar Foam::particle::trackToFace
622 (
623  const vector& displacement,
624  const scalar fraction
625 )
626 {
627  scalar f = 1;
628 
629  label tetTriI = onFace() ? 0 : -1;
630 
631  facei_ = -1;
632 
633  while (true)
634  {
635  f *= trackToTri(f*displacement, f*fraction, tetTriI);
636 
637  if (tetTriI == -1)
638  {
639  // The track has completed within the current tet
640  return 0;
641  }
642  else if (tetTriI == 0)
643  {
644  // The track has hit a face, so set the current face and return
645  facei_ = tetFacei_;
646  return f;
647  }
648  else
649  {
650  // Move to the next tet and continue the track
651  changeTet(tetTriI);
652  }
653  }
654 }
655 
656 
658 (
659  const vector& displacement,
660  const scalar fraction,
661  label& tetTriI
662 )
663 {
664  const vector x0 = position();
665  const vector x1 = displacement;
666  const barycentric y0 = coordinates_;
667 
668  if (debug)
669  {
670  Info<< "Particle " << origId() << endl
671  << "Tracking from " << x0 << " to " << x0 + x1 << endl;
672  }
673 
674  // Get the tet geometry
675  vector centre;
676  scalar detA;
678  stationaryTetReverseTransform(centre, detA, T);
679 
680  if (debug)
681  {
682  vector o, b, v1, v2;
683  stationaryTetGeometry(o, b, v1, v2);
684  Info<< "Tet points o=" << o << ", b=" << b
685  << ", v1=" << v1 << ", v2=" << v2 << endl
686  << "Tet determinant = " << detA << endl
687  << "Start local coordinates = " << y0 << endl;
688  }
689 
690  // Get the factor by which the displacement is increased
691  const scalar f = detA >= 0 ? 1 : negativeSpaceDisplacementFactor;
692 
693  // Calculate the local tracking displacement
694  barycentric Tx1(f*x1 & T);
695 
696  if (debug)
697  {
698  Info<< "Local displacement = " << Tx1 << "/" << detA << endl;
699  }
700 
701  // Calculate the hit fraction
702  label iH = -1;
703  scalar muH = std::isnormal(detA) && detA <= 0 ? VGREAT : 1/detA;
704  for (label i = 0; i < 4; ++ i)
705  {
706  if (std::isnormal(Tx1[i]) && Tx1[i] < 0)
707  {
708  scalar mu = - y0[i]/Tx1[i];
709 
710  if (debug)
711  {
712  Info<< "Hit on tet face " << i << " at local coordinate "
713  << y0 + mu*Tx1 << ", " << mu*detA*100 << "\% of the "
714  << "way along the track" << endl;
715  }
716 
717  if (0 <= mu && mu < muH)
718  {
719  iH = i;
720  muH = mu;
721  }
722  }
723  }
724 
725  // Set the new coordinates
726  barycentric yH = y0 + muH*Tx1;
727 
728  // Clamp to zero any negative coordinates generated by round-off error
729  for (label i = 0; i < 4; ++ i)
730  {
731  yH.replace(i, i == iH ? 0 : max(0, yH[i]));
732  }
733 
734  // Re-normalise if within the tet
735  if (iH == -1)
736  {
737  yH /= cmptSum(yH);
738  }
739 
740  // Set the new position and hit index
741  coordinates_ = yH;
742  tetTriI = iH;
743 
744  if (debug)
745  {
746  if (iH != -1)
747  {
748  Info<< "Track hit tet face " << iH << " first" << endl;
749  }
750  else
751  {
752  Info<< "Track hit no tet faces" << endl;
753  }
754  Info<< "End local coordinates = " << yH << endl
755  << "End global coordinates = " << position() << endl
756  << "Tracking displacement = " << position() - x0 << endl
757  << muH*detA*100 << "\% of the step from " << stepFraction_ << " to "
758  << stepFraction_ + fraction << " completed" << endl << endl;
759  }
760 
761  // Set the proportion of the track that has been completed
762  stepFraction_ += fraction*muH*detA;
763 
764  return iH != -1 ? 1 - muH*detA : 0;
765 }
766 
767 
769 (
770  const vector& displacement,
771  const scalar fraction,
772  label& tetTriI
773 )
774 {
775  const vector x0 = position();
776  const vector x1 = displacement;
777  const barycentric y0 = coordinates_;
778 
779  if (debug)
780  {
781  Info<< "Particle " << origId() << endl
782  << "Tracking from " << x0 << " to " << x0 + x1 << endl;
783  }
784 
785  // Get the tet geometry
786  Pair<vector> centre;
789  movingTetReverseTransform(fraction, centre, detA, T);
790 
791  // Get the factor by which the displacement is increased
792  const scalar f = detA[0] >= 0 ? 1 : negativeSpaceDisplacementFactor;
793 
794  // Get the relative global position
795  const vector x0Rel = x0 - centre[0];
796  const vector x1Rel = f*x1 - centre[1];
797 
798  // Form the determinant and hit equations
799  cubicEqn detAEqn(sqr(detA[0])*detA[3], detA[0]*detA[2], detA[1], 1);
800  const barycentric yC(1, 0, 0, 0);
801  const barycentric hitEqnA =
802  ((x1Rel & T[2]) + detA[3]*yC)*sqr(detA[0]);
803  const barycentric hitEqnB =
804  ((x1Rel & T[1]) + (x0Rel & T[2]) + detA[2]*yC)*detA[0];
805  const barycentric hitEqnC =
806  ((x1Rel & T[0]) + (x0Rel & T[1]) + detA[1]*yC);
807  const barycentric hitEqnD = y0;
808  FixedList<cubicEqn, 4> hitEqn;
809  forAll(hitEqn, i)
810  {
811  hitEqn[i] = cubicEqn(hitEqnA[i], hitEqnB[i], hitEqnC[i], hitEqnD[i]);
812  }
813 
814  // Calculate the hit fraction
815  label iH = -1;
816  scalar muH = std::isnormal(detA[0]) && detA[0] <= 0 ? VGREAT : 1/detA[0];
817  for (label i = 0; i < 4; ++ i)
818  {
819  const Roots<3> mu = hitEqn[i].roots();
820 
821  for (label j = 0; j < 3; ++ j)
822  {
823  if (mu.type(j) == roots::real && hitEqn[i].derivative(mu[j]) < 0)
824  {
825  if (0 <= mu[j] && mu[j] < muH)
826  {
827  iH = i;
828  muH = mu[j];
829  }
830  }
831  }
832  }
833 
834  // Set the new coordinates
835  barycentric yH
836  (
837  hitEqn[0].value(muH),
838  hitEqn[1].value(muH),
839  hitEqn[2].value(muH),
840  hitEqn[3].value(muH)
841  );
842  // !!! <-- This fails if the tet collapses onto the particle, as detA tends
843  // to zero at the hit. In this instance, we can differentiate the hit and
844  // detA polynomials to find a limiting location, but this will not be on a
845  // triangle. We will then need to do a second track through the degenerate
846  // tet to find the final hit position. This second track is over zero
847  // distance and therefore can be of the static mesh type. This has not yet
848  // been implemented.
849  const scalar detAH = detAEqn.value(muH);
850  if (!std::isnormal(detAH))
851  {
853  << "A moving tet collapsed onto a particle. This is not supported. "
854  << "The mesh is too poor, or the motion too severe, for particle "
855  << "tracking to function." << exit(FatalError);
856  }
857  yH /= detAH;
858 
859  // Clamp to zero any negative coordinates generated by round-off error
860  for (label i = 0; i < 4; ++ i)
861  {
862  yH.replace(i, i == iH ? 0 : max(0, yH[i]));
863  }
864 
865  // Re-normalise if within the tet
866  if (iH == -1)
867  {
868  yH /= cmptSum(yH);
869  }
870 
871  // Set the new position and hit index
872  coordinates_ = yH;
873  tetTriI = iH;
874 
875  // Set the proportion of the track that has been completed
876  stepFraction_ += fraction*muH*detA[0];
877 
878  if (debug)
879  {
880  if (iH != -1)
881  {
882  Info<< "Track hit tet face " << iH << " first" << endl;
883  }
884  else
885  {
886  Info<< "Track hit no tet faces" << endl;
887  }
888  Info<< "End local coordinates = " << yH << endl
889  << "End global coordinates = " << position() << endl;
890  }
891 
892  return iH != -1 ? 1 - muH*detA[0] : 0;
893 }
894 
895 
896 Foam::scalar Foam::particle::trackToTri
897 (
898  const vector& displacement,
899  const scalar fraction,
900  label& tetTriI
901 )
902 {
903  if (mesh_.moving())
904  {
905  return trackToMovingTri(displacement, fraction, tetTriI);
906  }
907  else
908  {
909  return trackToStationaryTri(displacement, fraction, tetTriI);
910  }
911 }
912 
913 
915 {
916  const Vector<label>& dirs = mesh_.geometricD();
917 
918  bool isConstrained = false;
919  forAll(dirs, dirI)
920  {
921  if (dirs[dirI] == -1)
922  {
923  isConstrained = true;
924  break;
925  }
926  }
927 
928  if (isConstrained)
929  {
930  vector pos = position();
932  track(pos - position(), 0);
933  }
934 }
935 
936 
938 {
939  if (!onBoundaryFace())
940  {
942  << "Patch data was requested for a particle that isn't on a patch"
943  << exit(FatalError);
944  }
945 
946  if (mesh_.moving())
947  {
948  Pair<vector> centre, base, vertex1, vertex2;
949  movingTetGeometry(1, centre, base, vertex1, vertex2);
950 
951  n = triPointRef(base[0], vertex1[0], vertex2[0]).normal();
952  n /= mag(n);
953 
954  // Interpolate the motion of the three face vertices to the current
955  // coordinates
956  U =
957  coordinates_.b()*base[1]
958  + coordinates_.c()*vertex1[1]
959  + coordinates_.d()*vertex2[1];
960 
961  // The movingTetGeometry method gives the motion as a displacement
962  // across the time-step, so we divide by the time-step to get velocity
963  U /= mesh_.time().deltaTValue();
964  }
965  else
966  {
967  vector centre, base, vertex1, vertex2;
968  stationaryTetGeometry(centre, base, vertex1, vertex2);
969 
970  n = triPointRef(base, vertex1, vertex2).normal();
971  n /= mag(n);
972 
973  U = Zero;
974  }
975 }
976 
977 
979 {}
980 
981 
983 {}
984 
985 
986 Foam::scalar Foam::particle::wallImpactDistance(const vector&) const
987 {
988  return 0.0;
989 }
990 
991 
993 (
994  const vectorTensorTransform& transform
995 )
996 {
997  // Get the transformed position
998  const vector pos = transform.invTransformPosition(position());
999 
1000  // Break the topology
1001  celli_ = -1;
1002  tetFacei_ = -1;
1003  tetPti_ = -1;
1004  facei_ = -1;
1005 
1006  // Store the position in the barycentric data
1007  coordinates_ = barycentric(1 - cmptSum(pos), pos.x(), pos.y(), pos.z());
1008 
1009  // Transform the properties
1010  transformProperties(- transform.t());
1011  if (transform.hasR())
1012  {
1013  transformProperties(transform.R().T());
1014  }
1015 }
1016 
1017 
1019 {
1020  // Get the position from the barycentric data
1021  const vector pos(coordinates_.b(), coordinates_.c(), coordinates_.d());
1022 
1023  // Create some arbitrary topology for the supplied cell
1024  celli_ = celli;
1025  tetFacei_ = mesh_.cells()[celli_][0];
1026  tetPti_ = 1;
1027  facei_ = -1;
1028 
1029  // Get the reverse transform and directly set the coordinates from the
1030  // position. This isn't likely to be correct; the particle is probably not
1031  // in this tet. It will, however, generate the correct vector when the
1032  // position method is called. A referred particle should never be tracked,
1033  // so this approximate topology is good enough. By using the nearby cell we
1034  // minimise the error associated with the incorrect topology.
1035  coordinates_ = barycentric(1, 0, 0, 0);
1036  if (mesh_.moving())
1037  {
1038  Pair<vector> centre;
1039  FixedList<scalar, 4> detA;
1041  movingTetReverseTransform(0, centre, detA, T);
1042  coordinates_ += (pos - centre[0]) & T[0]/detA[0];
1043  }
1044  else
1045  {
1046  vector centre;
1047  scalar detA;
1049  stationaryTetReverseTransform(centre, detA, T);
1050  coordinates_ += (pos - centre) & T/detA;
1051  }
1052 }
1053 
1054 
1057  const polyMesh& procMesh,
1058  const label procCell,
1059  const label procTetFace
1060 ) const
1061 {
1062  // The tet point on the procMesh differs from the current tet point if the
1063  // mesh and procMesh faces are of differing orientation. The change is the
1064  // same as in particle::correctAfterParallelTransfer.
1065 
1066  if
1067  (
1068  (mesh_.faceOwner()[tetFacei_] == celli_)
1069  == (procMesh.faceOwner()[procTetFace] == procCell)
1070  )
1071  {
1072  return tetPti_;
1073  }
1074  else
1075  {
1076  return procMesh.faces()[procTetFace].size() - 1 - tetPti_;
1077  }
1078 }
1079 
1080 
1083  const vector& position,
1084  const mapPolyMesh& mapper
1085 )
1086 {
1087  locate
1088  (
1089  position,
1090  nullptr,
1091  mapper.reverseCellMap()[celli_],
1092  true,
1093  "Particle mapped to a location outside of the mesh."
1094  );
1095 }
1096 
1097 
1098 // * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
1099 
1100 bool Foam::operator==(const particle& pA, const particle& pB)
1101 {
1102  return (pA.origProc() == pB.origProc() && pA.origId() == pB.origId());
1103 }
1104 
1105 
1106 bool Foam::operator!=(const particle& pA, const particle& pB)
1107 {
1108  return !(pA == pB);
1109 }
1110 
1111 
1112 // ************************************************************************* //
label origProc() const
Return the originating processor ID.
Definition: particleI.H:93
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
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
bool moving() const
Is mesh moving.
Definition: polyMesh.H:496
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
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:192
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 approproate for decomposition or reconstruction.
Definition: particle.C:1056
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 also stops on internal faces.
Definition: particle.C:622
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 also stops on tet triangles. On.
Definition: particle.C:897
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:418
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:509
vector invTransformPosition(const vector &v) const
Inverse transform the given position.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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:46
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:978
void replace(const direction, const Cmpt &)
Definition: VectorSpaceI.H:159
scalar f1
Definition: createFields.H:28
void prepareForInteractionListReferral(const vectorTensorTransform &transform)
Break the topology and store the particle position so that the.
Definition: particle.C:993
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:81
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:807
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:116
dimensionedScalar pos(const dimensionedScalar &ds)
Cubic equation of the form a*x^3 + b*x^2 + c*x + d = 0.
Definition: cubicEqn.H:49
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particle.C:937
void Swap(T &a, T &b)
Definition: Swap.H:43
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:524
3D tensor transformation operations.
face triFace(3)
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
Definition: particle.C:1018
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1049
const Cmpt & d() const
Definition: BarycentricI.H:80
static const zero Zero
Definition: zero.H:91
virtual scalar wallImpactDistance(const vector &n) const
The nearest distance to a wall that.
Definition: particle.C:986
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1036
const Cmpt & x() const
Definition: VectorI.H:75
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:105
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1082
const Cmpt & c() const
Definition: BarycentricI.H:73
bool onFace() const
Is the particle on a face?
Definition: particleI.H:180
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:526
const Time & time() const
Return time.
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:769
#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.
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
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:377
Field< vector > vectorField
Specialisation of Field<T> for vector.
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:186
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:658
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:603
bool operator!=(const particle &, const particle &)
Definition: particle.C:1106
void constrainToMeshCentre()
Set the constrained components of the particle position to the.
Definition: particle.C:914
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:204
Namespace for OpenFOAM.