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