particle.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 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 "tracking.H"
28 #include "polyTopoChangeMap.H"
29 #include "transform.H"
30 #include "treeDataCell.H"
31 #include "indexedOctree.H"
32 #include "cubicEqn.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
37 
38 namespace Foam
39 {
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const polyMesh& mesh,
49  const barycentric& coordinates,
50  const label celli,
51  const label tetFacei,
52  const label tetPti,
53  const label facei
54 )
55 :
56  coordinates_(coordinates),
57  celli_(celli),
58  tetFacei_(tetFacei),
59  tetPti_(tetPti),
60  facei_(facei),
61  stepFraction_(1),
62  stepFractionBehind_(0),
63  nTracksBehind_(0),
64  origProc_(Pstream::myProcNo()),
65  origId_(getNewParticleIndex())
66 {}
67 
68 
70 (
71  const polyMesh& mesh,
72  const vector& position,
73  const label celli,
74  label& nLocateBoundaryHits
75 )
76 :
77  coordinates_(- vGreat, - vGreat, - vGreat, - vGreat),
78  celli_(celli),
79  tetFacei_(-1),
80  tetPti_(-1),
81  facei_(-1),
82  stepFraction_(1),
83  stepFractionBehind_(0),
84  nTracksBehind_(0),
85  origProc_(Pstream::myProcNo()),
86  origId_(getNewParticleIndex())
87 {
88  if (!locate(mesh, position, celli))
89  {
90  nLocateBoundaryHits ++;
91  }
92 }
93 
94 
96 :
97  coordinates_(p.coordinates_),
98  celli_(p.celli_),
99  tetFacei_(p.tetFacei_),
100  tetPti_(p.tetPti_),
101  facei_(p.facei_),
102  stepFraction_(p.stepFraction_),
103  stepFractionBehind_(p.stepFractionBehind_),
104  nTracksBehind_(p.nTracksBehind_),
105  origProc_(p.origProc_),
106  origId_(p.origId_)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 (
114  const polyMesh& mesh,
115  const vector& position,
116  label celli
117 )
118 {
119  celli_ = celli;
120 
121  return
123  (
124  mesh, position,
125  coordinates_, celli_, tetFacei_, tetPti_, 1,
126  debug
127  ? static_cast<const string&>(string("Particle " + name(origId())))
128  : NullObjectRef<string>()
129  );
130 }
131 
132 
134 (
135  const polyMesh& mesh,
136  const vector& displacement,
137  const scalar fraction
138 )
139 {
140  const Tuple2<bool, scalar> onBoundaryAndF =
142  (
143  mesh, displacement, fraction,
144  coordinates_, celli_, tetFacei_, tetPti_, stepFraction_,
145  stepFractionBehind_, nTracksBehind_,
146  debug
147  ? static_cast<const string&>(string("Particle " + name(origId())))
148  : NullObjectRef<string>()
149  );
150 
151  facei_ = onBoundaryAndF.first() ? tetFacei_ : -1;
152 
153  return onBoundaryAndF.second();
154 }
155 
156 
158 (
159  const polyMesh& mesh,
160  const vector& displacement,
161  const scalar fraction
162 )
163 {
164  const Tuple2<bool, scalar> inNextCellAndF =
166  (
167  mesh, displacement, fraction,
168  coordinates_, celli_, tetFacei_, tetPti_, stepFraction_,
169  stepFractionBehind_, nTracksBehind_,
170  debug
171  ? static_cast<const string&>(string("Particle " + name(origId())))
172  : NullObjectRef<string>()
173  );
174 
175  facei_ = inNextCellAndF.first() ? tetFacei_ : -1;
176 
177  return inNextCellAndF.second();
178 }
179 
180 
182 (
183  const polyMesh& mesh,
184  const vector& displacement,
185  const scalar fraction
186 )
187 {
188  const Tuple2<bool, scalar> onFaceAndF =
190  (
191  mesh, displacement, fraction,
192  coordinates_, celli_, tetFacei_, tetPti_, stepFraction_,
193  stepFractionBehind_, nTracksBehind_,
194  debug
195  ? static_cast<const string&>(string("Particle " + name(origId())))
196  : NullObjectRef<string>()
197  );
198 
199  facei_ = onFaceAndF.first() ? tetFacei_ : -1;
200 
201  return onFaceAndF.second();
202 }
203 
204 
206 (
207  const polyMesh& mesh
208 ) const
209 {
210  if (cmptMin(mesh.geometricD()) == -1)
211  {
212  vector pos = position(mesh), posC = pos;
214  return pos - posC;
215  }
216  else
217  {
218  return vector::zero;
219  }
220 }
221 
222 
224 {}
225 
226 
228 {
229  const processorPolyPatch& ppp =
230  refCast<const processorPolyPatch>
231  (
233  );
234 
235  tracking::inProcessor(ppp, celli_, tetFacei_);
236 }
237 
238 
240 {
241  const processorPolyPatch& ppp =
242  refCast<const processorPolyPatch>
243  (
244  td.mesh.boundaryMesh()[td.sendToPatch]
245  );
246 
247  tracking::outProcessor(ppp, coordinates_, celli_, tetFacei_, tetPti_);
248 
249  // Remain on the face
250  facei_ = tetFacei_;
251 
252  // Transform the properties
253  if (ppp.transform().transformsPosition())
254  {
255  transformProperties(ppp.transform());
256  }
257 }
258 
259 
261 (
262  const polyMesh& mesh,
263  const label sendFromPatch,
264  const label sendToPatchFace,
265  const vector& sendToPosition
266 )
267 {
268  const nonConformalCyclicPolyPatch& nccpp =
269  static_cast<const nonConformalCyclicPolyPatch&>
270  (
271  mesh.boundaryMesh()[sendFromPatch]
272  );
273 
274  // Store the position in the barycentric data
275  coordinates_ =
277  (
278  1 - cmptSum(sendToPosition),
279  sendToPosition.x(),
280  sendToPosition.y(),
281  sendToPosition.z()
282  );
283 
284  // Break the topology
285  celli_ = -1;
286  tetFacei_ = -1;
287  tetPti_ = -1;
288 
289  // Store the local patch face in the face index
290  facei_ = sendToPatchFace;
291 
292  // Transform the properties
293  if (nccpp.transform().transformsPosition())
294  {
295  transformProperties(nccpp.nbrPatch().transform());
296  }
297 }
298 
299 
301 (
302  const polyMesh& mesh,
303  const label sendToPatch,
304  labelList& patchNLocateBoundaryHits
305 )
306 {
307  const nonConformalCyclicPolyPatch& nccpp =
308  static_cast<const nonConformalCyclicPolyPatch&>
309  (
310  mesh.boundaryMesh()[sendToPatch]
311  );
312 
313  // Get the position from the barycentric data
314  const vector receivePos
315  (
316  coordinates_.b(),
317  coordinates_.c(),
318  coordinates_.d()
319  );
320 
321  // Locate the particle on the receiving side
322  const label celli = mesh.faceOwner()[facei_ + nccpp.origPatch().start()];
323  if (!locate(mesh, receivePos, celli))
324  {
325  patchNLocateBoundaryHits[sendToPatch] ++;
326  }
327 
328  // The particle must remain associated with a face for the tracking to
329  // register as incomplete
330  facei_ = tetFacei_;
331 }
332 
333 
335 (
336  const polyMesh& mesh,
337  const transformer& transform
338 )
339 {
340  // Get the transformed position
341  const vector pos = transform.invTransformPosition(position(mesh));
342 
343  // Break the topology
344  celli_ = -1;
345  tetFacei_ = -1;
346  tetPti_ = -1;
347 
348  // Store the position in the barycentric data
349  coordinates_ = barycentric(1 - cmptSum(pos), pos.x(), pos.y(), pos.z());
350 
351  // Transform the properties
352  if (transform.transformsPosition())
353  {
354  transformProperties(inv(transform));
355  }
356 }
357 
358 
360 (
361  const polyMesh& mesh,
362  const label celli
363 )
364 {
365  // Get the position from the barycentric data
366  const vector pos(coordinates_.b(), coordinates_.c(), coordinates_.d());
367 
368  // Create some arbitrary topology for the supplied cell
369  celli_ = celli;
370  tetFacei_ = mesh.cells()[celli_][0];
371  tetPti_ = 1;
372 
373  // Calculate the coordinates for the assumed topology. This isn't likely to
374  // be correct; the particle is probably not in this tet. It will, however,
375  // generate the correct vector when the position method is called. A
376  // referred particle should never be tracked, so this approximate topology
377  // is good enough. By using the nearby cell we minimise the error
378  // associated with the incorrect topology.
379  coordinates_ =
381  (
382  mesh, pos,
383  celli_, tetFacei_, tetPti_, stepFraction_
384  );
385 }
386 
387 
389 (
390  const polyMesh& mesh,
391  const polyMesh& procMesh,
392  const label procCell,
393  const label procTetFace
394 ) const
395 {
396  // The tet point on the procMesh differs from the current tet point if the
397  // mesh and procMesh faces are of differing orientation. The change is the
398  // same as in particle::correctAfterParallelTransfer.
399 
400  if
401  (
402  (mesh.faceOwner()[tetFacei_] == celli_)
403  == (procMesh.faceOwner()[procTetFace] == procCell)
404  )
405  {
406  return tetPti_;
407  }
408  else
409  {
410  return procMesh.faces()[procTetFace].size() - 1 - tetPti_;
411  }
412 }
413 
414 
415 // * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
416 //
417 
418 bool Foam::operator==(const particle& pA, const particle& pB)
419 {
420  return (pA.origProc() == pB.origProc() && pA.origId() == pB.origId());
421 }
422 
423 
424 bool Foam::operator!=(const particle& pA, const particle& pB)
425 {
426  return !(pA == pB);
427 }
428 
429 
430 // ************************************************************************* //
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Inter-processor communications stream.
Definition: Pstream.H:56
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
static const Form zero
Definition: VectorSpace.H:118
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
virtual const transformer & transform() const
Inherit the cyclic transform method.
const nonConformalCyclicPolyPatch & nbrPatch() const
Neighbour patch.
const polyPatch & origPatch() const
Original patch.
label sendFromPatch
Patch from which to send the particle.
Definition: particle.H:112
label sendToPatch
Patch to which to send the particle.
Definition: particle.H:115
const polyMesh & mesh
Reference to the mesh.
Definition: particle.H:102
Base particle class.
Definition: particle.H:83
label origProc() const
Return the originating processor ID.
Definition: particleI.H:96
void prepareForNonConformalCyclicTransfer(const polyMesh &mesh, const label sendToPatch, const label sendToPatchFace, const vector &sendToPosition)
Make changes prior to a transfer across a non conformal cyclic.
Definition: particle.C:261
scalar trackToCell(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a new cell is reached.
Definition: particle.C:158
void prepareForProcessorTransfer(trackingData &td)
Make changes prior to a transfer across a processor boundary.
Definition: particle.C:227
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label facei)
Construct from components.
Definition: particle.C:47
void correctAfterInteractionListReferral(const polyMesh &mesh, const label celli)
Correct the topology after referral. Locates the particle.
Definition: particle.C:360
void correctAfterNonConformalCyclicTransfer(const polyMesh &mesh, const label sendToPatch, labelList &patchNLocateBoundaryHits)
Make changes following a transfer across a non conformal cyclic.
Definition: particle.C:301
scalar track(const polyMesh &mesh, const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:134
scalar trackToFace(const polyMesh &mesh, const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:182
void prepareForInteractionListReferral(const polyMesh &mesh, const transformer &transform)
Break the topology and store the cartesian position so that the.
Definition: particle.C:335
bool locate(const polyMesh &mesh, const vector &position, label celli)
Locate the particle at the given position.
Definition: particle.C:113
vector deviationFromMeshCentre(const polyMesh &mesh) const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:206
void correctAfterProcessorTransfer(trackingData &td)
Make changes following a transfer across a processor boundary.
Definition: particle.C:239
label procTetPt(const polyMesh &mesh, const polyMesh &procMesh, const label procCell, const label procTetFace) const
Return the tet point appropriate for decomposition or.
Definition: particle.C:389
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:159
virtual void transformProperties(const transformer &)
Transform the physical properties of the particle.
Definition: particle.C:223
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:108
static label particleCount
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:203
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1012
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
const cellList & cells() const
Neighbour processor patch.
virtual const transformer & transform() const
Return null transform between processor patches.
Vector-tensor class used to perform translations, rotations and scaling operations in 3D space.
Definition: transformer.H:84
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:613
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
Tuple2< bool, scalar > toCell(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
As toFace, except that if the track ends on an internal face then this.
Tuple2< bool, scalar > toBoundary(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
As toFace, except that the track continues across multiple cells until.
void outProcessor(const processorPolyPatch &outPatch, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Complete crossing of a processor patch. Restore the topology.
Definition: tracking.C:1777
void inProcessor(const processorPolyPatch &inPatch, label &celli, label &facei)
Initialise crossing of a processor patch. Breaks the topology in order.
Definition: tracking.C:1762
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
Tuple2< bool, scalar > toFace(const polyMesh &mesh, const Displacement &displacement, const scalar fraction, barycentric &coordinates, label &celli, label &facei, label &faceTrii, scalar &stepFraction, scalar &stepFractionBehind, label &nTracksBehind, const string &debugPrefix=NullObjectRef< string >())
Track along the displacement for a given fraction of the overall.
bool locate(const polyMesh &mesh, const point &position, barycentric &coordinates, label &celli, label &facei, label &faceTrii, const scalar stepFraction, const string &debugPrefix=NullObjectRef< string >())
Initialise the location at the given position. Returns whether or not a.
Definition: tracking.C:1593
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
bool operator!=(const particle &, const particle &)
Definition: particle.C:424
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:45
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:501
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:507
volScalarField & p
3D tensor transformation operations.