BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / lagrangian / basic / particle / particle.H
blob08913b41459c8d21c57210fc1157d318d42ec834
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Class
25     Foam::particle
27 Description
28     Base particle class
30 \*---------------------------------------------------------------------------*/
32 #ifndef particle_H
33 #define particle_H
35 #include "vector.H"
36 #include "Cloud.H"
37 #include "IDLList.H"
38 #include "labelList.H"
39 #include "pointField.H"
40 #include "faceList.H"
41 #include "OFstream.H"
42 #include "tetPointRef.H"
43 #include "FixedList.H"
44 #include "polyMeshTetDecomposition.H"
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 namespace Foam
51 // Forward declaration of classes
52 class particle;
54 class polyPatch;
56 class cyclicPolyPatch;
57 class processorPolyPatch;
58 class symmetryPolyPatch;
59 class wallPolyPatch;
60 class wedgePolyPatch;
62 // Forward declaration of friend functions and operators
64 Ostream& operator<<
66     Ostream&,
67     const particle&
70 bool operator==(const particle&, const particle&);
72 bool operator!=(const particle&, const particle&);
74 /*---------------------------------------------------------------------------*\
75                           Class Particle Declaration
76 \*---------------------------------------------------------------------------*/
78 class particle
80     public IDLList<particle>::link
82 public:
84     template<class CloudType>
85     class TrackingData
86     {
87         // Private data
89             //- Reference to the cloud containing (this) particle
90             CloudType& cloud_;
93     public:
95         // Public data
97             typedef CloudType cloudType;
99             //- Flag to switch processor
100             bool switchProcessor;
102             //- Flag to indicate whether to keep particle (false = delete)
103             bool keepParticle;
106         // Constructor
107         TrackingData(CloudType& cloud)
108         :
109             cloud_(cloud)
110         {}
113         // Member functions
115             //- Return a reference to the cloud
116             CloudType& cloud()
117             {
118                 return cloud_;
119             }
120     };
123 protected:
125     // Protected data
127         //- Reference to the polyMesh database
128         const polyMesh& mesh_;
130         //- Position of particle
131         vector position_;
133         //- Index of the cell it is in
134         label cellI_;
136         //- Face index if the particle is on a face otherwise -1
137         label faceI_;
139         //- Fraction of time-step completed
140         scalar stepFraction_;
142         //- Index of the face that owns the decomposed tet that the
143         //  particle is in
144         label tetFaceI_;
146         //- Index of the point on the face that defines the decomposed
147         //  tet that the particle is in.  Relative to the face base
148         //  point.
149         label tetPtI_;
151         //- Originating processor id
152         label origProc_;
154         //- Local particle id on originating processor
155         label origId_;
158     // Private Member Functions
160         //- Find the tet tri faces between position and tet centre
161         void findTris
162         (
163             const vector& position,
164             DynamicList<label>& faceList,
165             const tetPointRef& tet,
166             const FixedList<vector, 4>& tetAreas,
167             const FixedList<label, 4>& tetPlaneBasePtIs,
168             const scalar tol
169         ) const;
171         //- Find the lambda value for the line to-from across the
172         //  given tri face, where p = from + lambda*(to - from)
173         inline scalar tetLambda
174         (
175             const vector& from,
176             const vector& to,
177             const label triI,
178             const vector& tetArea,
179             const label tetPlaneBasePtI,
180             const label cellI,
181             const label tetFaceI,
182             const label tetPtI,
183             const scalar tol
184         ) const;
186         //- Find the lambda value for a moving tri face
187         inline scalar movingTetLambda
188         (
189             const vector& from,
190             const vector& to,
191             const label triI,
192             const vector& tetArea,
193             const label tetPlaneBasePtI,
194             const label cellI,
195             const label tetFaceI,
196             const label tetPtI,
197             const scalar tol
198         ) const;
200         //- Modify the tet owner data by crossing triI
201         inline void tetNeighbour(label triI);
203         //- Cross the from the given face across the given edge of the
204         //  given cell to find the resulting face and tetPtI
205         inline void crossEdgeConnectedFace
206         (
207             const label& cellI,
208             label& tetFaceI,
209             label& tetPtI,
210             const edge& e
211         );
213         //- Hit wall faces in the current cell if the
214         //- wallImpactDistance is non-zero.  They may not be in
215         //- different tets to the current.
216         template<class CloudType>
217         inline void hitWallFaces
218         (
219             const CloudType& td,
220             const vector& from,
221             const vector& to,
222             scalar& lambdaMin,
223             tetIndices& closestTetIs
224         );
227     // Patch interactions
229         //- Overridable function to handle the particle hitting a face
230         template<class TrackData>
231         void hitFace(TrackData& td);
233         //- Overridable function to handle the particle hitting a
234         //  patch.  Executed before other patch-hitting functions.
235         //  trackFraction is passed in to allow mesh motion to
236         //  interpolate in time to the correct face state.
237         template<class TrackData>
238         bool hitPatch
239         (
240             const polyPatch&,
241             TrackData& td,
242             const label patchI,
243             const scalar trackFraction,
244             const tetIndices& tetIs
245         );
247         //- Overridable function to handle the particle hitting a wedgePatch
248         template<class TrackData>
249         void hitWedgePatch(const wedgePolyPatch&, TrackData& td);
251         //- Overridable function to handle the particle hitting a
252         //  symmetryPatch
253         template<class TrackData>
254         void hitSymmetryPatch(const symmetryPolyPatch&, TrackData& td);
256         //- Overridable function to handle the particle hitting a cyclicPatch
257         template<class TrackData>
258         void hitCyclicPatch(const cyclicPolyPatch&, TrackData& td);
260         //- Overridable function to handle the particle hitting a
261         //  processorPatch
262         template<class TrackData>
263         void hitProcessorPatch(const processorPolyPatch&, TrackData& td);
265         //- Overridable function to handle the particle hitting a wallPatch
266         template<class TrackData>
267         void hitWallPatch
268         (
269             const wallPolyPatch&,
270             TrackData& td,
271             const tetIndices& tetIs
272         );
274         //- Overridable function to handle the particle hitting a
275         //  general patch
276         template<class TrackData>
277         void hitPatch(const polyPatch&, TrackData& td);
280 public:
282     // Static data members
284         //- Runtime type information
285         TypeName("particle");
287         //- String representation of properties
288         static string propHeader;
290         //- Cumulative particle counter - used to provode unique ID
291         static label particleCount_;
293         //- Fraction of distance to tet centre to move a particle to
294         // 'rescue' it from a tracking problem
295         static const scalar trackingCorrectionTol;
297         //- Fraction of the cell volume to use in determining tolerance values
298         //  for the denominator and numerator of lambda
299         static const scalar lambdaDistanceToleranceCoeff;
302     // Constructors
304         //- Construct from components
305         particle
306         (
307             const polyMesh& mesh,
308             const vector& position,
309             const label cellI,
310             const label tetFaceI,
311             const label tetPtI
312         );
314         //- Construct from components, tetFaceI_ and tetPtI_ are not
315         //  supplied so they will be deduced by a search
316         particle
317         (
318             const polyMesh& mesh,
319             const vector& position,
320             const label cellI,
321             bool doCellFacePt = true
322         );
324         //- Construct from Istream
325         particle(const polyMesh& mesh, Istream&, bool readFields = true);
327         //- Construct as a copy
328         particle(const particle& p);
330         //- Construct as a copy with refernce to a new mesh
331         particle(const particle& p, const polyMesh& mesh);
333         //- Construct a clone
334         virtual autoPtr<particle> clone() const
335         {
336             return autoPtr<particle>(new particle(*this));
337         }
340         //- Factory class to read-construct particles used for
341         //  parallel transfer
342         class iNew
343         {
344             const polyMesh& mesh_;
346         public:
348             iNew(const polyMesh& mesh)
349             :
350                 mesh_(mesh)
351             {}
353             autoPtr<particle> operator()(Istream& is) const
354             {
355                 return autoPtr<particle>(new particle(mesh_, is, true));
356             }
357         };
360     //- Destructor
361     virtual ~particle()
362     {}
365     // Member Functions
367         // Access
369             //- Get unique particle creation id
370             inline label getNewParticleID() const;
372             //- Return the mesh database
373             inline const polyMesh& mesh() const;
375             //- Return current particle position
376             inline const vector& position() const;
378             //- Return current particle position
379             inline vector& position();
381             //- Return current cell particle is in
382             inline label& cell();
384             //- Return current cell particle is in
385             inline label cell() const;
387             //- Return current tet face particle is in
388             inline label& tetFace();
390             //- Return current tet face particle is in
391             inline label tetFace() const;
393             //- Return current tet face particle is in
394             inline label& tetPt();
396             //- Return current tet face particle is in
397             inline label tetPt() const;
399             //- Return the indices of the current tet that the
400             //  particle occupies.
401             inline tetIndices currentTetIndices() const;
403             //- Return the geometry of the current tet that the
404             //  particle occupies.
405             inline tetPointRef currentTet() const;
407             //- Return the normal of the tri on tetFaceI_ for the
408             //  current tet.
409             inline vector normal() const;
411             //- Return the normal of the tri on tetFaceI_ for the
412             //  current tet at the start of the timestep, i.e. based
413             //  on oldPoints
414             inline vector oldNormal() const;
416             //- Return current face particle is on otherwise -1
417             inline label& face();
419             //- Return current face particle is on otherwise -1
420             inline label face() const;
422             //- Return the impact model to be used, soft or hard (default).
423             inline bool softImpact() const;
425             //- Return the particle current time
426             inline scalar currentTime() const;
429         // Check
431             //- Check the stored cell value (setting if necessary) and
432             //  initialise the tetFace and tetPt values
433             inline void initCellFacePt();
435             //- Is the particle on the boundary/(or outside the domain)?
436             inline bool onBoundary() const;
438             //- Is this global face an internal face?
439             inline bool internalFace(const label faceI) const;
441             //- Is this global face a boundary face?
442             inline bool boundaryFace(const label faceI) const;
444             //- Which patch is particle on
445             inline label patch(const label faceI) const;
447             //- Which face of this patch is this particle on
448             inline label patchFace
449             (
450                 const label patchI,
451                 const label faceI
452             ) const;
454             //- Return the fraction of time-step completed
455             inline scalar& stepFraction();
457             //-  Return the fraction of time-step completed
458             inline scalar stepFraction() const;
460             //- Return const access to the originating processor id
461             inline label origProc() const;
463             //- Return the originating processor id for manipulation
464             inline label& origProc();
466             //- Return const access to  the particle id on originating processor
467             inline label origId() const;
469             //- Return the particle id on originating processor for manipulation
470             inline label& origId();
473         // Track
475             //- Track particle to end of trajectory
476             //  or until it hits the boundary.
477             //  On entry 'stepFraction()' should be set to the fraction of the
478             //  time-step at which the tracking starts and on exit it contains
479             //  the fraction of the time-step completed.
480             //  Returns the boundary face index if the track stops at the
481             //  boundary, -1 otherwise.
482             template<class TrackData>
483             label track(const vector& endPosition, TrackData& td);
485             //- Track particle to a given position and returns 1.0 if the
486             //  trajectory is completed without hitting a face otherwise
487             //  stops at the face and returns the fraction of the trajectory
488             //  completed.
489             //  on entry 'stepFraction()' should be set to the fraction of the
490             //  time-step at which the tracking starts.
491             template<class TrackData>
492             scalar trackToFace(const vector& endPosition, TrackData& td);
494             //- Return the index of the face to be used in the interpolation
495             //  routine
496             inline label faceInterpolation() const;
499     // Transformations
501         //- Transform the physical properties of the particle
502         //  according to the given transformation tensor
503         virtual void transformProperties(const tensor& T);
505         //- Transform the physical properties of the particle
506         //  according to the given separation vector
507         virtual void transformProperties(const vector& separation);
509         //- The nearest distance to a wall that
510         //  the particle can be in the n direction
511         virtual scalar wallImpactDistance(const vector& n) const;
514     // Parallel transfer
516         //- Convert global addressing to the processor patch
517         //  local equivalents
518         template<class TrackData>
519         void prepareForParallelTransfer(const label patchI, TrackData& td);
521         //- Convert processor patch addressing to the global equivalents
522         //  and set the cellI to the face-neighbour
523         template<class TrackData>
524         void correctAfterParallelTransfer(const label patchI, TrackData& td);
527     // I-O
529         //- Read the fields associated with the owner cloud
530         template<class CloudType>
531         static void readFields(CloudType& c);
533         //- Write the fields associated with the owner cloud
534         template<class CloudType>
535         static void writeFields(const CloudType& c);
537         //- Write the particle data
538         void write(Ostream& os, bool writeFields) const;
541     // Friend Operators
543         friend Ostream& operator<<(Ostream&, const particle&);
545         friend bool operator==(const particle& pA, const particle& pB);
547         friend bool operator!=(const particle& pA, const particle& pB);
551 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
553 } // End namespace Foam
555 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
557 #include "particleI.H"
559 #define defineParticleTypeNameAndDebug(Type, DebugSwitch)                     \
560     template<>                                                                \
561     const Foam::word Particle<Type>::typeName(#Type);                         \
562     template<>                                                                \
563     int Particle<Type>::debug(Foam::debug::debugSwitch(#Type, DebugSwitch));
565 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
567 #ifdef NoRepository
568 #   include "particleTemplates.C"
569 #endif
571 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
573 #endif
575 // ************************************************************************* //