particle.H
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, 2020 OpenFOAM Foundation
9 Copyright (C) 2017-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::particle
29
30Description
31 Base particle class
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef particle_H
36#define particle_H
37
38#include "vector.H"
39#include "barycentric.H"
40#include "barycentricTensor.H"
41#include "Cloud.H"
42#include "IDLList.H"
43#include "pointField.H"
44#include "faceList.H"
45#include "OFstream.H"
46#include "tetPointRef.H"
47#include "FixedList.H"
49#include "particleMacros.H"
51
52// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57// Forward Declarations
58class particle;
59class polyPatch;
60class cyclicPolyPatch;
61class cyclicAMIPolyPatch;
62class cyclicACMIPolyPatch;
63class processorPolyPatch;
64class symmetryPlanePolyPatch;
65class symmetryPolyPatch;
66class wallPolyPatch;
67class wedgePolyPatch;
68
69Ostream& operator<<(Ostream&, const particle&);
70bool operator==(const particle&, const particle&);
71bool operator!=(const particle&, const particle&);
72
73/*---------------------------------------------------------------------------*\
74 Class Particle Declaration
75\*---------------------------------------------------------------------------*/
77class particle
78:
79 public IDLList<particle>::link
80{
81 // Private Data
82
83 //- Size in bytes of the position data
84 static const std::size_t sizeofPosition;
85
86 //- Size in bytes of the fields
87 static const std::size_t sizeofFields;
88
89 //- The value of nBehind_ at which tracking is abandoned. See the
90 // description of nBehind_.
91 static const label maxNBehind_;
92
93
94public:
96 class trackingData
97 {
98 public:
99
100 // Public Data
101
102 //- Flag to switch processor
103 bool switchProcessor;
104
105 //- Flag to indicate whether to keep particle (false = delete)
106 bool keepParticle;
107
108
109 // Constructor
110 template <class TrackCloudType>
111 trackingData(const TrackCloudType& cloud)
112 {}
113 };
114
115
116 //- Old particle positions content for OpenFOAM-1706 and earlier
120 label celli;
121 label facei;
123 label tetFacei;
124 label tetPti;
125 label origProc;
126 label origId;
127 };
128
129
130private:
131
132 // Private Data
133
134 //- Reference to the polyMesh database
135 const polyMesh& mesh_;
136
137 //- Coordinates of particle
138 barycentric coordinates_;
139
140 //- Index of the cell it is in
141 label celli_;
142
143 //- Index of the face that owns the decomposed tet that the
144 //- particle is in
145 label tetFacei_;
146
147 //- Index of the point on the face that defines the decomposed
148 //- tet that the particle is in. Relative to the face base
149 //- point.
150 label tetPti_;
151
152 //- Face index if the particle is on a face otherwise -1
153 label facei_;
154
155 //- Fraction of time-step completed
156 scalar stepFraction_;
157
158 //- The distance behind the maximum distance reached so far
159 scalar behind_;
160
161 //- The number of tracks carried out that ended in a distance behind the
162 //- maximum distance reached so far. Once this reaches maxNBehind_,
163 // tracking is abandoned for the current step. This is needed because
164 // when tetrahedra are inverted a straight trajectory can form a closed
165 // loop through regions of overlapping positive and negative space.
166 // Without this break clause, such loops can result in a valid track
167 // which never ends.
168 label nBehind_;
169
170 //- Originating processor id
171 label origProc_;
172
173 //- Local particle id on originating processor
174 label origId_;
175
176
177 // Private Member Functions
178
179 // Tetrahedra functions
180
181 //- Get the vertices of the current tet
182 inline void stationaryTetGeometry
183 (
184 vector& centre,
185 vector& base,
186 vector& vertex1,
187 vector& vertex2
188 ) const;
189
190 //- Get the transformation associated with the current tet. This
191 // will convert a barycentric position within the tet to a
192 // Cartesian position in the global coordinate system. The
193 // conversion is x = A & y, where x is the Cartesian position, y is
194 // the barycentric position and A is the transformation tensor.
195 inline barycentricTensor stationaryTetTransform() const;
196
197 //- Get the reverse transform associated with the current tet. The
198 // conversion is detA*y = (x - centre) & T. The variables x, y and
199 // centre have the same meaning as for the forward transform. T is
200 // the transposed inverse of the forward transform tensor, A,
201 // multiplied by its determinant, detA. This separation allows
202 // the barycentric tracking algorithm to function on inverted or
203 // degenerate tetrahedra.
204 void stationaryTetReverseTransform
205 (
206 vector& centre,
207 scalar& detA,
209 ) const;
210
211 //- Get the vertices of the current moving tet. Two values are
212 // returned for each vertex. The first is a constant, and the
213 // second is a linear coefficient of the track fraction.
214 inline void movingTetGeometry
215 (
216 const scalar endStepFraction,
217 Pair<vector>& centre,
218 Pair<vector>& base,
219 Pair<vector>& vertex1,
220 Pair<vector>& vertex2
221 ) const;
222
223 //- Get the transformation associated with the current, moving, tet.
224 // This is of the same form as for the static case. As with the
225 // moving geometry, a linear function of the tracking fraction is
226 // returned for each component.
227 inline Pair<barycentricTensor> movingTetTransform
228 (
229 const scalar endStepFraction
230 ) const;
231
232 //- Get the reverse transformation associated with the current,
233 // moving, tet. This is of the same form as for the static case. As
234 // with the moving geometry, a function of the tracking fraction is
235 // returned for each component. The functions are higher order than
236 // for the forward transform; the determinant is cubic, and the
237 // tensor is quadratic.
238 void movingTetReverseTransform
239 (
240 const scalar endStepFraction,
241 Pair<vector>& centre,
244 ) const;
245
246
247 // Transformations
248
249 //- Reflection transform. Corrects the coordinates when the particle
250 // moves between two tets which share a base vertex, but for which
251 // the other two non cell-centre vertices are reversed. All hits
252 // which retain the same face behave this way, as do face hits.
253 void reflect();
254
255 //- Rotation transform. Corrects the coordinates when the particle
256 // moves between two tets with different base vertices, but are
257 // otherwise similarly oriented. Hits which change the face within
258 // the cell make use of both this and the reflect transform.
259 void rotate(const bool direction);
260
261
262 // Topology changes
263
264 //- Change tet within a cell. Called after a triangle is hit.
265 void changeTet(const label tetTriI);
266
267 //- Change tet face within a cell. Called by changeTet.
268 void changeFace(const label tetTriI);
269
270 //- Change cell. Called when the particle hits an internal face.
271 void changeCell();
272
273 //- Put the particle on the lowest indexed patch for the current set
274 // of coincident faces. In the case of an ACMI-wall pair, this will
275 // move the particle from the wall face to the ACMI face, because
276 // ACMI patches are always listed before their associated non-
277 // overlapping patch.
278 void changeToMasterPatch();
279
280
281 // Geometry changes
282
283 //- Locate the particle at the given position
284 void locate
285 (
286 const vector& position,
287 const vector* direction,
288 const label celli,
289 const bool boundaryFail,
290 const string& boundaryMsg
291 );
292
293
294protected:
295
296 // Patch interactions
297
298 //- Read particle from stream. Optionally (for old format) return
299 // read position. Used by construct-from-Istream
300 void readData
301 (
302 Istream& is,
304 const bool readFields,
305 const bool newFormat,
306 const bool doLocate
307 );
308
309 //- Overridable function to handle the particle hitting a patch.
310 // Executed before other patch-hitting functions.
311 template<class TrackCloudType>
312 bool hitPatch(TrackCloudType&, trackingData&);
313
314 //- Overridable function to handle the particle hitting a wedgePatch
315 template<class TrackCloudType>
316 void hitWedgePatch(TrackCloudType&, trackingData&);
317
318 //- Overridable function to handle the particle hitting a
319 // symmetryPlanePatch
320 template<class TrackCloudType>
321 void hitSymmetryPlanePatch(TrackCloudType&, trackingData&);
322
323 //- Overridable function to handle the particle hitting a symmetryPatch
324 template<class TrackCloudType>
325 void hitSymmetryPatch(TrackCloudType&, trackingData&);
326
327 //- Overridable function to handle the particle hitting a cyclicPatch
328 template<class TrackCloudType>
329 void hitCyclicPatch(TrackCloudType&, trackingData&);
330
331 //- Overridable function to handle the particle hitting a cyclicAMIPatch
332 template<class TrackCloudType>
333 void hitCyclicAMIPatch(TrackCloudType&, trackingData&, const vector&);
334
335 //- Overridable function to handle the particle hitting a
336 // cyclicACMIPatch
337 template<class TrackCloudType>
338 void hitCyclicACMIPatch(TrackCloudType&, trackingData&, const vector&);
339
340 //- Overridable function to handle the particle hitting a processorPatch
341 template<class TrackCloudType>
342 void hitProcessorPatch(TrackCloudType&, trackingData&);
343
344 //- Overridable function to handle the particle hitting a wallPatch
345 template<class TrackCloudType>
346 void hitWallPatch(TrackCloudType&, trackingData&);
347
348
349 //- Dispatch function for boundary face interaction. Calls one of
350 // the above (hitWedgePatch, hitCyclicPatch etc) depending on the
351 // patch type
352 template<class TrackCloudType>
353 void hitBoundaryFace
354 (
355 const vector& direction,
356 TrackCloudType& cloud,
357 trackingData& td
358 );
359
360
361public:
362
363 // Static Data Members
364
365 //- Runtime type information
366 TypeName("particle");
367
368 //- String representation of properties
370 (
371 "(coordinatesa coordinatesb coordinatesc coordinatesd) "
372 "celli tetFacei tetPti facei stepFraction origProc origId"
373 );
374
375 //- Cumulative particle counter - used to provide unique ID
376 static label particleCount_;
377
378 //- Write particle coordinates file (v1712 and later)
379 //- Default is true
380 static bool writeLagrangianCoordinates;
381
382 //- Write particle positions file (v1706 format and earlier)
383 //- Default is true (disable in etc/controlDict)
384 static bool writeLagrangianPositions;
385
386
387 // Constructors
388
389 //- Construct from components
391 (
392 const polyMesh& mesh,
394 const label celli,
395 const label tetFacei,
396 const label tetPti
397 );
398
399 //- Construct from a position and a cell.
400 // Searches for the rest of the required topology.
402 (
403 const polyMesh& mesh,
404 const vector& position,
405 const label celli = -1
406 );
407
408 //- Construct from position components
410 (
411 const polyMesh& mesh,
412 const vector& position,
413 const label celli,
414 const label tetFacei,
415 const label tetPti,
416 const bool doLocate = true
417 );
418
419
420 //- Construct from Istream
422 (
423 const polyMesh& mesh,
424 Istream&,
425 const bool readFields = true,
426 const bool newFormat = true,
427 const bool doLocate = true
428 );
429
430 //- Construct as a copy
431 particle(const particle& p);
432
433 //- Construct as a copy with reference to a new mesh
434 particle(const particle& p, const polyMesh& mesh);
435
436 //- Construct a clone
437 virtual autoPtr<particle> clone() const
438 {
439 return autoPtr<particle>::New(*this);
440 }
441
442 //- Factory class to read-construct particles (for parallel transfer)
443 class iNew
444 {
445 const polyMesh& mesh_;
446
447 public:
449 iNew(const polyMesh& mesh)
450 :
451 mesh_(mesh)
452 {}
455 {
456 return autoPtr<particle>::New(mesh_, is, true);
457 }
458 };
459
460
461 //- Destructor
462 virtual ~particle() = default;
463
464
465 // Member Functions
466
467 // Access
468
469 //- Get unique particle creation id
470 inline label getNewParticleID() const;
471
472 //- Return the mesh database
473 inline const polyMesh& mesh() const noexcept;
474
475 //- Return current particle coordinates
476 inline const barycentric& coordinates() const noexcept;
477
478 //- Return current cell particle is in
479 inline label cell() const noexcept;
480
481 //- Return current cell particle is in for manipulation
482 inline label& cell() noexcept;
483
484 //- Return current tet face particle is in
485 inline label tetFace() const noexcept;
486
487 //- Return current tet face particle is in for manipulation
488 inline label& tetFace() noexcept;
489
490 //- Return current tet face particle is in
491 inline label tetPt() const noexcept;
492
493 //- Return current tet face particle is in for manipulation
494 inline label& tetPt() noexcept;
495
496 //- Return current face particle is on otherwise -1
497 inline label face() const noexcept;
498
499 //- Return current face particle is on for manipulation
500 inline label& face() noexcept;
501
502 //- Return the fraction of time-step completed
503 inline scalar stepFraction() const noexcept;
504
505 //- Return the fraction of time-step completed
506 inline scalar& stepFraction() noexcept;
507
508 //- Return the originating processor ID
509 inline label origProc() const noexcept;
510
511 //- Return the originating processor ID
512 inline label& origProc() noexcept;
513
514 //- Return the particle ID on the originating processor
515 inline label origId() const noexcept;
516
517 //- Return the particle ID on the originating processor
518 inline label& origId() noexcept;
519
520
521 // Check
522
523 //- Return the step fraction change within the overall time-step.
524 // Returns the start value and the change as a scalar pair. Always
525 // return Pair<scalar>(0, 1), unless sub-cycling is in effect, in
526 // which case the values will reflect the span of the sub-cycle
527 // within the time-step.
528 inline Pair<scalar> stepFractionSpan() const;
529
530 //- Return the current fraction within the timestep. This differs
531 // from the stored step fraction due to sub-cycling.
532 inline scalar currentTimeFraction() const;
533
534 //- Return indices of the current tet that the particle occupies.
536
537 //- Return the current tet transformation tensor
539
540 //- The (unit) normal of the tri on tetFacei_ for the current tet.
541 inline vector normal() const;
542
543 //- Is the particle on a face?
544 inline bool onFace() const noexcept;
545
546 //- Is the particle on an internal face?
547 inline bool onInternalFace() const noexcept;
548
549 //- Is the particle on a boundary face?
550 inline bool onBoundaryFace() const noexcept;
551
552 //- Return the index of patch that the particle is on
553 inline label patch() const;
554
555 //- Return current particle position
556 inline vector position() const;
557
558 //- Reset particle data
559 inline void reset();
560
561
562 // Track
563
564 //- Track along the displacement for a given fraction of the overall
565 // step. End when the track is complete, or when a boundary is hit.
566 // On exit, stepFraction_ will have been incremented to the current
567 // position, and facei_ will be set to the index of the boundary
568 // face that was hit, or -1 if the track completed within a cell.
569 // The proportion of the displacement still to be completed is
570 // returned.
571 scalar track
572 (
573 const vector& displacement,
574 const scalar fraction
575 );
576
577 //- As particle::track, but also stops on internal faces.
578 scalar trackToFace
579 (
580 const vector& displacement,
581 const scalar fraction
582 );
583
584 //- As particle::trackToFace, but also stops on tet triangles. On
585 // exit, tetTriI is set to the index of the tet triangle that was
586 // hit, or -1 if the end position was reached.
587 scalar trackToTri
588 (
589 const vector& displacement,
590 const scalar fraction,
591 label& tetTriI
592 );
593
594 //- As particle::trackToTri, but for stationary meshes
596 (
597 const vector& displacement,
598 const scalar fraction,
599 label& tetTriI
600 );
601
602 //- As particle::trackToTri, but for moving meshes
603 scalar trackToMovingTri
604 (
605 const vector& displacement,
606 const scalar fraction,
607 label& tetTriI
608 );
609
610 //- Hit the current face. If the current face is internal than this
611 // crosses into the next cell. If it is a boundary face then this will
612 // interact the particle with the relevant patch.
613 template<class TrackCloudType>
614 void hitFace
615 (
616 const vector& direction,
617 TrackCloudType& cloud,
618 trackingData& td
619 );
620
621 //- Convenience function. Combines trackToFace and hitFace
622 template<class TrackCloudType>
624 (
625 const vector& direction,
626 const scalar fraction,
627 TrackCloudType& cloud,
628 trackingData& td
629 );
630
631 //- Get the displacement from the mesh centre. Used to correct the
632 // particle position in cases with reduced dimensionality. Returns a
633 // zero vector for three-dimensional cases.
635
636
637 // Patch data
638
639 //- Get the normal and velocity of the current patch location
640 inline void patchData(vector& n, vector& U) const;
641
642
643 // Transformations
644
645 //- Transform the physical properties of the particle
646 // according to the given transformation tensor
647 virtual void transformProperties(const tensor& T);
648
649 //- Transform the physical properties of the particle
650 // according to the given separation vector
651 virtual void transformProperties(const vector& separation);
652
653
654 // Parallel transfer
655
656 //- Convert global addressing to the processor patch local equivalents
658
659 //- Convert processor patch addressing to the global equivalents
660 // and set the celli to the face-neighbour
661 void correctAfterParallelTransfer(const label patchi, trackingData& td);
662
663
664 // Interaction list referral
665
666 //- Break the topology and store the particle position so that the
667 // particle can be referred.
669 (
671 );
672
673 //- Correct the topology after referral. The particle may still be
674 // outside the stored tet and therefore not track-able.
675 void correctAfterInteractionListReferral(const label celli);
676
677
678 // Decompose and reconstruct
679
680 //- Return the tet point appropriate for decomposition or reconstruction
681 // to or from the given mesh.
682 label procTetPt
683 (
684 const polyMesh& procMesh,
685 const label procCell,
686 const label procTetFace
687 ) const;
688
689
690 // Mapping
691
692 //- Map after a topology change
693 void autoMap(const vector& position, const mapPolyMesh& mapper);
694
695 //- Set the addressing based on the provided position
696 void relocate(const point& position, const label celli = -1);
697
698
699
700 // I-O
701
702 //- Write the name representation to stream
703 template<class Type>
704 static void writePropertyName
705 (
706 Ostream& os,
707 const word& name,
708 const word& delim
709 );
710
711 //- Write a named particle property to stream,
712 //- optionally filtered based on its name
713 template<class Type>
714 static void writeProperty
715 (
716 Ostream& os,
717 const word& name,
718 const Type& value,
719 const bool nameOnly,
720 const word& delim,
721 const wordRes& filters = wordRes::null()
722 );
723
724 //- Write a named particle property list to stream,
725 //- optionally filtered based on its name
726 template<class Type>
727 static void writeProperty
728 (
729 Ostream& os,
730 const word& name,
731 const Field<Type>& values,
732 const bool nameOnly,
733 const word& delim,
734 const wordRes& filters = wordRes::null()
735 );
736
737 //- Read the fields associated with the owner cloud
738 template<class TrackCloudType>
739 static void readFields(TrackCloudType& c);
740
741 //- Write the fields associated with the owner cloud
742 template<class TrackCloudType>
743 static void writeFields(const TrackCloudType& c);
744
745 //- Write individual particle properties to stream
746 void writeProperties
747 (
748 Ostream& os,
749 const wordRes& filters,
750 const word& delim,
751 const bool namesOnly
752 ) const;
753
754 //- Read particle fields as objects from the obr registry
755 template<class CloudType>
756 static void readObjects(CloudType& c, const objectRegistry& obr);
757
758 //- Write particle fields as objects into the obr registry
759 // Always writes "position", not "coordinate"
760 template<class CloudType>
761 static void writeObjects(const CloudType& c, objectRegistry& obr);
762
763 //- Write the particle barycentric coordinates and cell info
764 void writeCoordinates(Ostream& os) const;
765
766 //- Write the particle position and cell id
767 virtual void writePosition(Ostream& os) const;
768
769
770 // Friend Operators
772 friend Ostream& operator<<(Ostream&, const particle&);
774 friend bool operator==(const particle& pA, const particle& pB);
776 friend bool operator!=(const particle& pA, const particle& pB);
777};
778
779
780// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
781
782} // End namespace Foam
783
784// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
785
786#include "particleI.H"
787
788// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
789
790#ifdef NoRepository
791 #include "particleTemplates.C"
792#endif
793
794// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
795
796#endif
797
798// ************************************************************************* //
Intrusive doubly-linked list.
label n
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Generic templated field type.
Definition: Field.H:82
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
Template class for intrusive linked lists.
Definition: ILList.H:69
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:69
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
Factory class to read-construct particles (for parallel transfer)
Definition: particle.H:443
iNew(const polyMesh &mesh)
Definition: particle.H:448
autoPtr< particle > operator()(Istream &is) const
Definition: particle.H:453
bool switchProcessor
Flag to switch processor.
Definition: particle.H:102
trackingData(const TrackCloudType &cloud)
Definition: particle.H:110
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
Base particle class.
Definition: particle.H:79
void correctAfterInteractionListReferral(const label celli)
Correct the topology after referral. The particle may still be.
Definition: particle.C:1162
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:1072
label procTetPt(const polyMesh &procMesh, const label procCell, const label procTetFace) const
Return the tet point appropriate for decomposition or reconstruction.
Definition: particle.C:1200
tetIndices currentTetIndices() const noexcept
Return indices of the current tet that the particle occupies.
Definition: particleI.H:265
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:1057
static void writeProperty(Ostream &os, const word &name, const Type &value, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:271
scalar trackToStationaryTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for stationary meshes.
Definition: particle.C:714
void correctAfterParallelTransfer(const label patchi, trackingData &td)
Convert processor patch addressing to the global equivalents.
Definition: particle.C:1088
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but also stops on internal faces.
Definition: particle.C:658
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a.
vector position() const
Return current particle position.
Definition: particleI.H:314
label tetPt() const noexcept
Return current tet face particle is in.
Definition: particleI.H:173
void hitFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
static void writeObjects(const CloudType &c, objectRegistry &obr)
Write particle fields as objects into the obr registry.
void prepareForInteractionListReferral(const vectorTensorTransform &transform)
Break the topology and store the particle position so that the.
Definition: particle.C:1137
void hitBoundaryFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Dispatch function for boundary face interaction. Calls one of.
bool onFace() const noexcept
Is the particle on a face?
Definition: particleI.H:290
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:137
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
label tetFace() const noexcept
Return current tet face particle is in.
Definition: particleI.H:161
void autoMap(const vector &position, const mapPolyMesh &mapper)
Map after a topology change.
Definition: particle.C:1226
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition: particleI.H:143
scalar trackToMovingTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToTri, but for moving meshes.
Definition: particle.C:846
virtual autoPtr< particle > clone() const
Construct a clone.
Definition: particle.H:436
bool onBoundaryFace() const noexcept
Is the particle on a boundary face?
Definition: particleI.H:302
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:308
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition: particleI.H:197
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
void readData(Istream &is, point &position, const bool readFields, const bool newFormat, const bool doLocate)
Read particle from stream. Optionally (for old format) return.
Definition: particleIO.C:76
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particleI.H:328
void prepareForParallelTransfer()
Convert global addressing to the processor patch local equivalents.
Definition: particle.C:1080
void hitCyclicAMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a cyclicAMIPatch.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
Pair< scalar > stepFractionSpan() const
Return the step fraction change within the overall time-step.
Definition: particleI.H:233
static label particleCount_
Cumulative particle counter - used to provide unique ID.
Definition: particle.H:375
static bool writeLagrangianPositions
Definition: particle.H:383
void writeProperties(Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly) const
Write individual particle properties to stream.
Definition: particleIO.C:220
scalar track(const vector &displacement, const scalar fraction)
Track along the displacement for a given fraction of the overall.
Definition: particle.C:639
TypeName("particle")
Runtime type information.
bool onInternalFace() const noexcept
Is the particle on an internal face?
Definition: particleI.H:296
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:123
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
label origId() const noexcept
Return the particle ID on the originating processor.
Definition: particleI.H:221
void reset()
Reset particle data.
Definition: particleI.H:320
vector normal() const
The (unit) normal of the tri on tetFacei_ for the current tet.
Definition: particleI.H:284
virtual void writePosition(Ostream &os) const
Write the particle position and cell id.
Definition: particleIO.C:264
void relocate(const point &position, const label celli=-1)
Set the addressing based on the provided position.
Definition: particle.C:1242
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.
Definition: particleI.H:257
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
scalar trackToTri(const vector &displacement, const scalar fraction, label &tetTriI)
As particle::trackToFace, but also stops on tet triangles. On.
Definition: particle.C:1040
static void writePropertyName(Ostream &os, const word &name, const word &delim)
Write the name representation to stream.
static bool writeLagrangianCoordinates
Definition: particle.H:379
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
label origProc() const noexcept
Return the originating processor ID.
Definition: particleI.H:209
void writeCoordinates(Ostream &os) const
Write the particle barycentric coordinates and cell info.
Definition: particleIO.C:245
virtual ~particle()=default
Destructor.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
Vector-tensor class used to perform translations and rotations in 3D space.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & T
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
uint8_t direction
Definition: direction.H:56
const direction noexcept
Definition: Scalar.H:223
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Macros for adding to particle property lists.
#define DefinePropertyList(str)
Define a static 'propertyList' for particle properties.
Old particle positions content for OpenFOAM-1706 and earlier.
Definition: particle.H:117
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73