fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / lagrangian / basic / Particle / Particle.H
blob643c6eb5d5060511fbdfa1d71dff511632a83a62
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::Particle
28 Description
30 \*---------------------------------------------------------------------------*/
32 #ifndef Particle_H
33 #define Particle_H
35 #include "vector.H"
36 #include "IDLList.H"
37 #include "labelList.H"
38 #include "pointField.H"
39 #include "faceList.H"
40 #include "typeInfo.H"
41 #include "OFstream.H"
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 namespace Foam
48 template<class Particle>
49 class Cloud;
51 class wedgePolyPatch;
52 class symmetryPolyPatch;
53 class cyclicPolyPatch;
54 class processorPolyPatch;
55 class wallPolyPatch;
56 class polyPatch;
58 // Forward declaration of friend functions and operators
60 template<class ParticleType>
61 class Particle;
63 template<class ParticleType>
64 Ostream& operator<<
66     Ostream&,
67     const Particle<ParticleType>&
71 /*---------------------------------------------------------------------------*\
72                            Class Particle Declaration
73 \*---------------------------------------------------------------------------*/
75 template<class ParticleType>
76 class Particle
78     public IDLList<ParticleType>::link
81 public:
83     //- Class used to pass tracking data to the trackToFace function
84     class trackData
85     {
87         // Private data
89             //- Reference to the cloud containing this particle
90             Cloud<ParticleType>& cloud_;
93     public:
95         bool switchProcessor;
96         bool keepParticle;
99         // Constructors
101             inline trackData(Cloud<ParticleType>& cloud);
104         // Member functions
106             //- Return a reference to the cloud
107             inline Cloud<ParticleType>& cloud();
108     };
111 protected:
113     // Private data
115         //- Reference to the particle cloud
116         const Cloud<ParticleType>& cloud_;
118         //- Position of particle
119         vector position_;
121         //- Index of the cell it is in
122         label celli_;
124         //- Face index if the particle is on a face otherwise -1
125         label facei_;
127         //- Fraction of time-step completed
128         scalar stepFraction_;
130         //- Originating processor id
131         label origProc_;
133         //- Local particle id on originating processor
134         label origId_;
137     // Private member functions
139         //- Return the 'lambda' value for the position, p, on the face,
140         // where, p = from + lamda*(to - from)
141         // for non-static meshes
142         inline scalar lambda
143         (
144             const vector& from,
145             const vector& to,
146             const label facei,
147             const scalar stepFraction
148         ) const;
150         //- Return the 'lambda' value for the position, p, on the face,
151         // where, p = from + lamda*(to - from)
152         // for static meshes
153         inline scalar lambda
154         (
155             const vector& from,
156             const vector& to,
157             const label facei
158         ) const;
160         //- Find the faces between position and cell centre
161         void findFaces
162         (
163             const vector& position,
164             DynamicList<label>& faceList
165         ) const;
167         //- Find the faces between position and cell centre
168         void findFaces
169         (
170             const vector& position,
171             const label celli,
172             const scalar stepFraction,
173             DynamicList<label>& faceList
174         ) const;
177     // Patch interactions
179         //- Overridable function to handle the particle hitting a patch
180         //  Executed before other patch-hitting functions
181         template<class TrackData>
182         bool hitPatch
183         (
184             const polyPatch&,
185             TrackData& td,
186             const label patchI
187         );
189         //- Overridable function to handle the particle hitting a wedgePatch
190         template<class TrackData>
191         void hitWedgePatch
192         (
193             const wedgePolyPatch&,
194             TrackData& td
195         );
197         //- Overridable function to handle the particle hitting a
198         //  symmetryPatch
199         template<class TrackData>
200         void hitSymmetryPatch
201         (
202             const symmetryPolyPatch&,
203             TrackData& td
204         );
206         //- Overridable function to handle the particle hitting a cyclicPatch
207         template<class TrackData>
208         void hitCyclicPatch
209         (
210             const cyclicPolyPatch&,
211             TrackData& td
212         );
214         //- Overridable function to handle the particle hitting a
215         //  processorPatch
216         template<class TrackData>
217         void hitProcessorPatch
218         (
219             const processorPolyPatch&,
220             TrackData& td
221         );
223         //- Overridable function to handle the particle hitting a wallPatch
224         template<class TrackData>
225         void hitWallPatch
226         (
227             const wallPolyPatch&,
228             TrackData& td
229         );
231         //- Overridable function to handle the particle hitting a
232         //  general patch
233         template<class TrackData>
234         void hitPatch
235         (
236             const polyPatch&,
237             TrackData& td
238         );
241     // Transformations
243         //- Transform the position the particle
244         //  according to the given transformation tensor
245         virtual void transformPosition(const tensor& T);
247         //- Transform the physical properties of the particle
248         //  according to the given transformation tensor
249         virtual void transformProperties(const tensor& T);
251         //- Transform the physical properties of the particle
252         //  according to the given separation vector
253         virtual void transformProperties(const vector& separation);
256     // Parallel transfer
258         //- Convert global addressing to the processor patch
259         //  local equivalents
260         template<class TrackData>
261         void prepareForParallelTransfer(const label patchi, TrackData& td);
263         //- Convert processor patch addressing to the global equivalents
264         //  and set the celli to the face-neighbour
265         template<class TrackData>
266         void correctAfterParallelTransfer(const label patchi, TrackData& td);
269 public:
271     friend class Cloud<ParticleType>;
274     // Static data members
276         //- Runtime type information
277         TypeName("Particle");
279         //- String representation of properties
280         static string propHeader;
283     // Constructors
285         //- Construct from components
286         Particle
287         (
288             const Cloud<ParticleType>&,
289             const vector& position,
290             const label celli
291         );
293         //- Construct from Istream
294         Particle
295         (
296             const Cloud<ParticleType>&,
297             Istream&,
298             bool readFields = true
299         );
301         //- Construct as a copy
302         Particle(const Particle& p);
304         //- Construct a clone
305         autoPtr<ParticleType> clone() const
306         {
307             return autoPtr<Particle>(new Particle(*this));
308         }
311         //- Factory class to read-construct particles used for
312         //  parallel transfer
313         class iNew
314         {
315             // Private data
317                 //- Reference to the cloud
318                 const Cloud<ParticleType>& cloud_;
321         public:
323             iNew(const Cloud<ParticleType>& cloud)
324             :
325                 cloud_(cloud)
326             {}
328             autoPtr<ParticleType> operator()(Istream& is) const
329             {
330                 return autoPtr<ParticleType>
331                 (
332                     new ParticleType
333                     (
334                         cloud_,
335                         is,
336                         true
337                     )
338                 );
339             }
340         };
343     //- Destructor
344     virtual ~Particle()
345     {}
348     // Member Functions
350         // Access
352             //- Return true if particle is in cell
353             inline bool inCell() const;
355             //- Return true if position is in cell i
356             inline bool inCell
357             (
358                 const vector& position,
359                 const label celli,
360                 const scalar stepFraction
361             ) const;
363             //- Return current particle position
364             inline const vector& position() const;
366             //- Return current particle position
367             inline vector& position();
369             //- Return current cell particle is in
370             inline label& cell();
372             //- Return current cell particle is in
373             inline label cell() const;
375             //- Return current face particle is on otherwise -1
376             inline label face() const;
378             //- Return reference to the particle cloud
379             inline const Cloud<ParticleType>& cloud() const;
381             //- Return the impact model to be used, soft or hard (default).
382             inline bool softImpact() const;
384             //- Return the particle current time
385             inline scalar currentTime() const;
388         // Check
390             //- Is the particle on the boundary/(or outside the domain)?
391             inline bool onBoundary() const;
393             //- Which patch is particle on
394             inline label patch(const label facei) const;
396             //- Which face of this patch is this particle on
397             inline label patchFace
398             (
399                 const label patchi,
400                 const label facei
401             ) const;
403             //- The nearest distance to a wall that
404             //  the particle can be in the n direction
405             inline scalar wallImpactDistance(const vector& n) const;
407             //- Return the fraction of time-step completed
408             inline scalar& stepFraction();
410             //-  Return the fraction of time-step completed
411             inline scalar stepFraction() const;
413             //- Return the originating processor id
414             inline label origProc() const;
416             //- Return the particle id on originating processor
417             inline label origId() const;
420         // Track
422             //- Track particle to end of trajectory
423             //  or until it hits the boundary.
424             //  On entry 'stepFraction()' should be set to the fraction of the
425             //  time-step at which the tracking starts and on exit it contains
426             //  the fraction of the time-step completed.
427             //  Returns the boundary face index if the track stops at the
428             //  boundary, -1 otherwise.
429             template<class TrackData>
430             label track
431             (
432                 const vector& endPosition,
433                 TrackData& td
434             );
436             //- Calls the templated track with dummy TrackData
437             label track(const vector& endPosition);
439             //- Track particle to a given position and returns 1.0 if the
440             //  trajectory is completed without hitting a face otherwise
441             //  stops at the face and returns the fraction of the trajectory
442             //  completed.
443             //  on entry 'stepFraction()' should be set to the fraction of the
444             //  time-step at which the tracking starts.
445             template<class TrackData>
446             scalar trackToFace
447             (
448                 const vector& endPosition,
449                 TrackData& td
450             );
452             //- Calls the templated trackToFace with dummy TrackData
453             scalar trackToFace(const vector& endPosition);
455             //- Return the index of the face to be used in the interpolation
456             //  routine
457             inline label faceInterpolation() const;
460     // I-O
462         //- Read the fields associated with the owner cloud
463         static void readFields(Cloud<ParticleType>& c);
465         //- Write the fields associated with the owner cloud
466         static void writeFields(const Cloud<ParticleType>& c);
468         //- Write the particle data
469         void write(Ostream& os, bool writeFields) const;
471     // Ostream Operator
473         friend Ostream& operator<< <ParticleType>
474         (
475             Ostream&,
476             const Particle<ParticleType>&
477         );
481 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
483 } // End namespace Foam
485 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
487 #include "ParticleI.H"
489 #define defineParticleTypeNameAndDebug(Type, DebugSwitch)                     \
490     template<>                                                                \
491     const Foam::word Particle<Type>::typeName(#Type);                         \
492     template<>                                                                \
493     int Particle<Type>::debug(Foam::debug::debugSwitch(#Type, DebugSwitch));
495 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
497 #ifdef NoRepository
498 #   include "Particle.C"
499 #endif
501 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
503 #endif
505 // ************************************************************************* //