| 12
 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
 505
 506
 507
 508
 509
 510
 511
 512
 513
 514
 515
 516
 517
 518
 519
 520
 521
 522
 523
 524
 525
 526
 527
 528
 529
 530
 531
 532
 533
 534
 535
 536
 537
 538
 539
 540
 541
 542
 543
 544
 545
 
 |  
#ifndef particle_H
#define particle_H
 
#include "vector.H"
#include "Cloud.H"
#include "IDLList.H"
#include "labelList.H"
#include "pointField.H"
#include "faceList.H"
#include "OFstream.H"
#include "tetPointRef.H"
#include "FixedList.H"
#include "polyMeshTetDecomposition.H"
 
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
namespace Foam
{
 
// Forward declaration of classes
class particle;
 
class polyPatch;
 
class cyclicPolyPatch;
class processorPolyPatch;
class symmetryPolyPatch;
class wallPolyPatch;
class wedgePolyPatch;
 
// Forward declaration of friend functions and operators
 
Ostream& operator<<
(
    Ostream&,
    const particle&
);
 
bool operator==(const particle&, const particle&);
 
bool operator!=(const particle&, const particle&);
 
/*---------------------------------------------------------------------------*\
                          Class Particle Declaration
\*---------------------------------------------------------------------------*/
 
class particle
:
    public IDLList<particle>::link
{
public:
 
    template<class CloudType>
    class TrackingData
    {
        // Private data
 
            //- Reference to the cloud containing (this) particle
            CloudType& cloud_;
 
 
    public:
 
        // Public data
 
            typedef CloudType cloudType;
 
            //- Flag to switch processor
            bool switchProcessor;
 
            //- Flag to indicate whether to keep particle (false = delete)
            bool keepParticle;
 
 
        // Constructor
        TrackingData(CloudType& cloud)
        :
            cloud_(cloud)
        {}
 
 
        // Member functions
 
            //- Return a reference to the cloud
            CloudType& cloud()
            {
                return cloud_;
            }
    };
 
 
protected:
 
    // Protected data
 
        //- Reference to the polyMesh database
        const polyMesh& mesh_;
 
        //- 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_;
 
        //- Index of the face that owns the decomposed tet that the
        //  particle is in
        label tetFaceI_;
 
        //- Index of the point on the face that defines the decomposed
        //  tet that the particle is in.  Relative to the face base
        //  point.
        label tetPtI_;
 
        //- Originating processor id
        label origProc_;
 
        //- Local particle id on originating processor
        label origId_;
 
 
    // Private Member Functions
 
        //- Find the tet tri faces between position and tet centre
        void findTris
        (
            const vector& position,
            DynamicList<label>& faceList,
            const tetPointRef& tet,
            const FixedList<vector, 4>& tetAreas,
            const FixedList<label, 4>& tetPlaneBasePtIs,
            const scalar tol
        ) const;
 
        //- Find the lambda value for the line to-from across the
        //  given tri face, where p = from + lambda*(to - from)
        inline scalar tetLambda
        (
            const vector& from,
            const vector& to,
            const label triI,
            const vector& tetArea,
            const label tetPlaneBasePtI,
            const label cellI,
            const label tetFaceI,
            const label tetPtI,
            const scalar tol
        ) const;
 
        //- Find the lambda value for a moving tri face
        inline scalar movingTetLambda
        (
            const vector& from,
            const vector& to,
            const label triI,
            const vector& tetArea,
            const label tetPlaneBasePtI,
            const label cellI,
            const label tetFaceI,
            const label tetPtI,
            const scalar tol
        ) const;
 
        //- Modify the tet owner data by crossing triI
        inline void tetNeighbour(label triI);
 
        //- Cross the from the given face across the given edge of the
        //  given cell to find the resulting face and tetPtI
        inline void crossEdgeConnectedFace
        (
            const label& cellI,
            label& tetFaceI,
            label& tetPtI,
            const edge& e
        );
 
        //- Hit wall faces in the current cell if the
        //- wallImpactDistance is non-zero.  They may not be in
        //- different tets to the current.
        template<class CloudType>
        inline void hitWallFaces
        (
            const CloudType& td,
            const vector& from,
            const vector& to,
            scalar& lambdaMin,
            tetIndices& closestTetIs
        );
 
 
    // Patch interactions
 
        //- Overridable function to handle the particle hitting a face
        template<class TrackData>
        void hitFace(TrackData& td);
 
        //- Overridable function to handle the particle hitting a
        //  patch.  Executed before other patch-hitting functions.
        //  trackFraction is passed in to allow mesh motion to
        //  interpolate in time to the correct face state.
        template<class TrackData>
        bool hitPatch
        (
            const polyPatch&,
            TrackData& td,
            const label patchI,
            const scalar trackFraction,
            const tetIndices& tetIs
        );
 
        //- 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,
            const tetIndices& tetIs
        );
 
        //- Overridable function to handle the particle hitting a
        //  general patch
        template<class TrackData>
        void hitPatch(const polyPatch&, TrackData& td);
 
 
public:
 
    // Static data members
 
        //- Runtime type information
        TypeName("particle");
 
        //- String representation of properties
        static string propHeader;
 
        //- Cumulative particle counter - used to provode unique ID
        static label particleCount_;
 
        //- Fraction of distance to tet centre to move a particle to
        // 'rescue' it from a tracking problem
        static const scalar trackingCorrectionTol;
 
        //- Fraction of the cell volume to use in determining tolerance values
        //  for the denominator and numerator of lambda
        static const scalar lambdaDistanceToleranceCoeff;
 
 
    // Constructors
 
        //- Construct from components
        particle
        (
            const polyMesh& mesh,
            const vector& position,
            const label cellI,
            const label tetFaceI,
            const label tetPtI
        );
 
        //- Construct from components, tetFaceI_ and tetPtI_ are not
        //  supplied so they will be deduced by a search
        particle
        (
            const polyMesh& mesh,
            const vector& position,
            const label cellI,
            bool doCellFacePt = true
        );
 
        //- Construct from Istream
        particle(const polyMesh& mesh, Istream&, bool readFields = true);
 
        //- Construct as a copy
        particle(const particle& p);
 
        //- Construct as a copy with refernce to a new mesh
        particle(const particle& p, const polyMesh& mesh);
 
        //- Construct a clone
        virtual autoPtr<particle> clone() const
        {
            return autoPtr<particle>(new particle(*this));
        }
 
 
        //- Factory class to read-construct particles used for
        //  parallel transfer
        class iNew
        {
            const polyMesh& mesh_;
 
        public:
 
            iNew(const polyMesh& mesh)
            :
                mesh_(mesh)
            {}
 
            autoPtr<particle> operator()(Istream& is) const
            {
                return autoPtr<particle>(new particle(mesh_, is, true));
            }
        };
 
 
    //- Destructor
    virtual ~particle()
    {}
 
 
    // Member Functions
 
        // Access
 
            //- Get unique particle creation id
            inline label getNewParticleID() const;
 
            //- Return the mesh database
            inline const polyMesh& mesh() 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 tet face particle is in
            inline label& tetFace();
 
            //- Return current tet face particle is in
            inline label tetFace() const;
 
            //- Return current tet face particle is in
            inline label& tetPt();
 
            //- Return current tet face particle is in
            inline label tetPt() const;
 
            //- Return the indices of the current tet that the
            //  particle occupies.
            inline tetIndices currentTetIndices() const;
 
            //- Return the geometry of the current tet that the
            //  particle occupies.
            inline tetPointRef currentTet() const;
 
            //- Return the normal of the tri on tetFaceI_ for the
            //  current tet.
            inline vector normal() const;
 
            //- Return the normal of the tri on tetFaceI_ for the
            //  current tet at the start of the timestep, i.e. based
            //  on oldPoints
            inline vector oldNormal() const;
 
            //- Return current face particle is on otherwise -1
            inline label& face();
 
            //- Return current face particle is on otherwise -1
            inline label face() 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
 
            //- Check the stored cell value (setting if necessary) and
            //  initialise the tetFace and tetPt values
            inline void initCellFacePt();
 
            //- Is the particle on the boundary/(or outside the domain)?
            inline bool onBoundary() const;
 
            //- Is this global face an internal face?
            inline bool internalFace(const label faceI) const;
 
            //- Is this global face a boundary face?
            inline bool boundaryFace(const label faceI) 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;
 
            //- Return the fraction of time-step completed
            inline scalar& stepFraction();
 
            //-  Return the fraction of time-step completed
            inline scalar stepFraction() const;
 
            //- Return const access to the originating processor id
            inline label origProc() const;
 
            //- Return the originating processor id for manipulation
            inline label& origProc();
 
            //- Return const access to  the particle id on originating processor
            inline label origId() const;
 
            //- Return the particle id on originating processor for manipulation
            inline label& origId();
 
 
        // 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);
 
            //- 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);
 
            //- Return the index of the face to be used in the interpolation
            //  routine
            inline label faceInterpolation() const;
 
 
    // Transformations
 
        //- 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);
 
        //- The nearest distance to a wall that
        //  the particle can be in the n direction
        virtual scalar wallImpactDistance(const vector& n) const;
 
 
    // 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);
 
 
    // I-O
 
        //- Read the fields associated with the owner cloud
        template<class CloudType>
        static void readFields(CloudType& c);
 
        //- Write the fields associated with the owner cloud
        template<class CloudType>
        static void writeFields(const CloudType& c);
 
        //- Write the particle data
        void write(Ostream& os, bool writeFields) const;
 
 
    // Friend Operators
 
        friend Ostream& operator<<(Ostream&, const particle&);
 
        friend bool operator==(const particle& pA, const particle& pB);
 
        friend bool operator!=(const particle& pA, const particle& pB);
};
 
 
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
} // End namespace Foam
 
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
#include "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 "particleTemplates.C"
#endif
 
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
 
#endif
 
// ************************************************************************* // | 
Partager