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