This file is indexed.

/usr/include/freefoam/lagrangian/Particle.H is in libfreefoam-dev 0.1.0+dfsg-1build1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Class
    Foam::Particle

Description

\*---------------------------------------------------------------------------*/

#ifndef Particle_H
#define Particle_H

#include <OpenFOAM/vector.H>
#include <OpenFOAM/IDLList.H>
#include <OpenFOAM/labelList.H>
#include <OpenFOAM/pointField.H>
#include <OpenFOAM/faceList.H>
#include <OpenFOAM/typeInfo.H>
#include <OpenFOAM/OFstream.H>

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

namespace Foam
{

template<class Particle>
class Cloud;

class wedgePolyPatch;
class symmetryPolyPatch;
class cyclicPolyPatch;
class processorPolyPatch;
class wallPolyPatch;
class polyPatch;

// Forward declaration of friend functions and operators

template<class ParticleType>
class Particle;

template<class ParticleType>
Ostream& operator<<
(
    Ostream&,
    const Particle<ParticleType>&
);


/*---------------------------------------------------------------------------*\
                           Class Particle Declaration
\*---------------------------------------------------------------------------*/

template<class ParticleType>
class Particle
:
    public IDLList<ParticleType>::link
{

public:

    //- Class used to pass tracking data to the trackToFace function
    class trackData
    {

        // Private data

            //- Reference to the cloud containing this particle
            Cloud<ParticleType>& cloud_;


    public:

        bool switchProcessor;
        bool keepParticle;


        // Constructors

            inline trackData(Cloud<ParticleType>& cloud);


        // Member functions

            //- Return a reference to the cloud
            inline Cloud<ParticleType>& cloud();
    };


protected:

    // Private data

        //- Reference to the particle cloud
        const Cloud<ParticleType>& cloud_;

        //- Position of particle
        vector position_;

        //- Index of the cell it is in
        label celli_;

        //- Face index if the particle is on a face otherwise -1
        label facei_;

        //- Fraction of time-step completed
        scalar stepFraction_;

        //- Originating processor id
        label origProc_;

        //- Local particle id on originating processor
        label origId_;


    // Private member functions

        //- Return the 'lambda' value for the position, p, on the face,
        // where, p = from + lamda*(to - from)
        // for non-static meshes
        inline scalar lambda
        (
            const vector& from,
            const vector& to,
            const label facei,
            const scalar stepFraction
        ) const;

        //- Return the 'lambda' value for the position, p, on the face,
        // where, p = from + lamda*(to - from)
        // for static meshes
        inline scalar lambda
        (
            const vector& from,
            const vector& to,
            const label facei
        ) const;

        //- Find the faces between position and cell centre
        void findFaces
        (
            const vector& position,
            DynamicList<label>& faceList
        ) const;

        //- Find the faces between position and cell centre
        void findFaces
        (
            const vector& position,
            const label celli,
            const scalar stepFraction,
            DynamicList<label>& faceList
        ) const;


    // Patch interactions

        //- Overridable function to handle the particle hitting a patch
        //  Executed before other patch-hitting functions
        template<class TrackData>
        bool hitPatch
        (
            const polyPatch&,
            TrackData& td,
            const label patchI
        );

        //- Overridable function to handle the particle hitting a wedgePatch
        template<class TrackData>
        void hitWedgePatch
        (
            const wedgePolyPatch&,
            TrackData& td
        );

        //- Overridable function to handle the particle hitting a
        //  symmetryPatch
        template<class TrackData>
        void hitSymmetryPatch
        (
            const symmetryPolyPatch&,
            TrackData& td
        );

        //- Overridable function to handle the particle hitting a cyclicPatch
        template<class TrackData>
        void hitCyclicPatch
        (
            const cyclicPolyPatch&,
            TrackData& td
        );

        //- Overridable function to handle the particle hitting a
        //  processorPatch
        template<class TrackData>
        void hitProcessorPatch
        (
            const processorPolyPatch&,
            TrackData& td
        );

        //- Overridable function to handle the particle hitting a wallPatch
        template<class TrackData>
        void hitWallPatch
        (
            const wallPolyPatch&,
            TrackData& td
        );

        //- Overridable function to handle the particle hitting a
        //  general patch
        template<class TrackData>
        void hitPatch
        (
            const polyPatch&,
            TrackData& td
        );


    // Transformations

        //- Transform the position the particle
        //  according to the given transformation tensor
        virtual void transformPosition(const tensor& T);

        //- Transform the physical properties of the particle
        //  according to the given transformation tensor
        virtual void transformProperties(const tensor& T);

        //- Transform the physical properties of the particle
        //  according to the given separation vector
        virtual void transformProperties(const vector& separation);


    // Parallel transfer

        //- Convert global addressing to the processor patch
        //  local equivalents
        template<class TrackData>
        void prepareForParallelTransfer(const label patchi, TrackData& td);

        //- Convert processor patch addressing to the global equivalents
        //  and set the celli to the face-neighbour
        template<class TrackData>
        void correctAfterParallelTransfer(const label patchi, TrackData& td);


public:

    friend class Cloud<ParticleType>;


    // Static data members

        //- Runtime type information
        TypeName("Particle");

        //- String representation of properties
        static string propHeader;


    // Constructors

        //- Construct from components
        Particle
        (
            const Cloud<ParticleType>&,
            const vector& position,
            const label celli
        );

        //- Construct from Istream
        Particle
        (
            const Cloud<ParticleType>&,
            Istream&,
            bool readFields = true
        );

        //- Construct as a copy
        Particle(const Particle& p);

        //- Construct a clone
        autoPtr<ParticleType> clone() const
        {
            return autoPtr<Particle>(new Particle(*this));
        }


        //- Factory class to read-construct particles used for
        //  parallel transfer
        class iNew
        {
            // Private data

                //- Reference to the cloud
                const Cloud<ParticleType>& cloud_;


        public:

            iNew(const Cloud<ParticleType>& cloud)
            :
                cloud_(cloud)
            {}

            autoPtr<ParticleType> operator()(Istream& is) const
            {
                return autoPtr<ParticleType>
                (
                    new ParticleType
                    (
                        cloud_,
                        is,
                        true
                    )
                );
            }
        };


    //- Destructor
    virtual ~Particle()
    {}


    // Member Functions

        // Access

            //- Return true if particle is in cell
            inline bool inCell() const;

            //- Return true if position is in cell i
            inline bool inCell
            (
                const vector& position,
                const label celli,
                const scalar stepFraction
            ) const;

            //- Return current particle position
            inline const vector& position() const;

            //- Return current particle position
            inline vector& position();

            //- Return current cell particle is in
            inline label& cell();

            //- Return current cell particle is in
            inline label cell() const;

            //- Return current face particle is on otherwise -1
            inline label face() const;

            //- Return reference to the particle cloud
            inline const Cloud<ParticleType>& cloud() const;

            //- Return the impact model to be used, soft or hard (default).
            inline bool softImpact() const;

            //- Return the particle current time
            inline scalar currentTime() const;


        // Check

            //- Is the particle on the boundary/(or outside the domain)?
            inline bool onBoundary() const;

            //- Which patch is particle on
            inline label patch(const label facei) const;

            //- Which face of this patch is this particle on
            inline label patchFace
            (
                const label patchi,
                const label facei
            ) const;

            //- The nearest distance to a wall that
            //  the particle can be in the n direction
            inline scalar wallImpactDistance(const vector& n) const;

            //- Return the fraction of time-step completed
            inline scalar& stepFraction();

            //-  Return the fraction of time-step completed
            inline scalar stepFraction() const;

            //- Return the originating processor id
            inline label origProc() const;

            //- Return the particle id on originating processor
            inline label origId() const;


        // Track

            //- Track particle to end of trajectory
            //  or until it hits the boundary.
            //  On entry 'stepFraction()' should be set to the fraction of the
            //  time-step at which the tracking starts and on exit it contains
            //  the fraction of the time-step completed.
            //  Returns the boundary face index if the track stops at the
            //  boundary, -1 otherwise.
            template<class TrackData>
            label track
            (
                const vector& endPosition,
                TrackData& td
            );

            //- Calls the templated track with dummy TrackData
            label track(const vector& endPosition);

            //- Track particle to a given position and returns 1.0 if the
            //  trajectory is completed without hitting a face otherwise
            //  stops at the face and returns the fraction of the trajectory
            //  completed.
            //  on entry 'stepFraction()' should be set to the fraction of the
            //  time-step at which the tracking starts.
            template<class TrackData>
            scalar trackToFace
            (
                const vector& endPosition,
                TrackData& td
            );

            //- Calls the templated trackToFace with dummy TrackData
            scalar trackToFace(const vector& endPosition);

            //- Return the index of the face to be used in the interpolation
            //  routine
            inline label faceInterpolation() const;


    // I-O

        //- Read the fields associated with the owner cloud
        static void readFields(Cloud<ParticleType>& c);

        //- Write the fields associated with the owner cloud
        static void writeFields(const Cloud<ParticleType>& c);

        //- Write the particle data
        void write(Ostream& os, bool writeFields) const;

    // Ostream Operator

        friend Ostream& operator<< <ParticleType>
        (
            Ostream&,
            const Particle<ParticleType>&
        );
};


// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

} // End namespace Foam

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#include <lagrangian/ParticleI.H>

#define defineParticleTypeNameAndDebug(Type, DebugSwitch)                     \
    template<>                                                                \
    const Foam::word Particle<Type>::typeName(#Type);                         \
    template<>                                                                \
    int Particle<Type>::debug(Foam::debug::debugSwitch(#Type, DebugSwitch));

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#ifdef NoRepository
#   include <lagrangian/Particle.C>
#endif

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

#endif

// ************************ vim: set sw=4 sts=4 et: ************************ //