ParticleI.H 31.7 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
6
7
8
9
10
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11
12
13
14
    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.
15
16
17
18
19
20
21

    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
22
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
23
24
25
26
27
28
29
30

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

#include "polyMesh.H"

// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

template<class ParticleType>
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
inline void Foam::Particle<ParticleType>::findTris
(
    const vector& position,
    DynamicList<label>& faceList,
    const tetPointRef& tet,
    const FixedList<vector, 4>& tetAreas,
    const FixedList<label, 4>& tetPlaneBasePtIs
) const
{
    faceList.clear();

    const point Ct = tet.centre();

    for (label i = 0; i < 4; i++)
    {
        scalar lambda = tetLambda
        (
            Ct,
            position,
            i,
            tetAreas[i],
            tetPlaneBasePtIs[i],
            cellI_,
            tetFaceI_,
            tetPtI_
        );

        if ((lambda > 0.0) && (lambda < 1.0))
        {
            faceList.append(i);
        }
    }
}


template<class ParticleType>
inline Foam::scalar Foam::Particle<ParticleType>::tetLambda
(
    const vector& from,
    const vector& to,
    const label triI,
    const vector& n,
    const label tetPlaneBasePtI,
    const label cellI,
    const label tetFaceI,
    const label tetPtI
) const
{
    const polyMesh& mesh = cloud_.polyMesh_;
    const pointField& pPts = mesh.points();

    if (mesh.moving())
    {
        return movingTetLambda
        (
            from,
            to,
            triI,
            n,
            tetPlaneBasePtI,
            cellI,
            tetFaceI,
            tetPtI
        );
    }

    const point& base = pPts[tetPlaneBasePtI];

    scalar lambdaNumerator = (base - from) & n;
    scalar lambdaDenominator = (to - from) & n;

    if (mag(lambdaDenominator) < SMALL)
    {
        if (mag(lambdaNumerator) < SMALL)
        {
            // Track starts on the face, and is potentially
            // parallel to it.  +-SMALL)/+-SMALL is not a good
            // comparison, return 0.0, in anticipation of tet
            // centre correction.

            return 0.0;
        }
        else
        {
            if (mag((to - from)) < SMALL)
            {
                // Zero length track, not along the face, face
                // cannot be crossed.

                return GREAT;
            }
            else
            {
                // Trajectory is non-zero and parallel to face

                lambdaDenominator = sign(lambdaDenominator)*SMALL;
            }
        }
    }

    return lambdaNumerator/lambdaDenominator;
}



template<class ParticleType>
inline Foam::scalar Foam::Particle<ParticleType>::movingTetLambda
138
139
140
(
    const vector& from,
    const vector& to,
141
142
143
144
145
146
    const label triI,
    const vector& n,
    const label tetPlaneBasePtI,
    const label cellI,
    const label tetFaceI,
    const label tetPtI
147
) const
148
149
{
    const polyMesh& mesh = cloud_.polyMesh_;
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
    const pointField& pPts = mesh.points();
    const pointField& oldPPts = mesh.oldPoints();

    // Base point of plane at end of motion
    const point& b = pPts[tetPlaneBasePtI];

    // n: Normal of plane at end of motion

    // Base point of plane at start of timestep
    const point& b00 = oldPPts[tetPlaneBasePtI];

    // Base point of plane at start of tracking portion (cast forward by
    // stepFraction)
    point b0 = b00 + stepFraction_*(b - b00);

    // Normal of plane at start of tracking portion
    vector n0 = vector::zero;
167
168

    {
169
170
171
172
        tetIndices tetIs(cellI, tetFaceI, tetPtI, mesh);

        // tet at timestep start
        tetPointRef tet00 = tetIs.oldTet(mesh);
173

174
175
176
177
178
179
180
181
182
183
184
185
        // tet at timestep end
        tetPointRef tet = tetIs.tet(mesh);

        point tet0PtA = tet00.a() + stepFraction_*(tet.a() - tet00.a());
        point tet0PtB = tet00.b() + stepFraction_*(tet.b() - tet00.b());
        point tet0PtC = tet00.c() + stepFraction_*(tet.c() - tet00.c());
        point tet0PtD = tet00.d() + stepFraction_*(tet.d() - tet00.d());

        // Tracking portion start tet (cast forward by stepFraction)
        tetPointRef tet0(tet0PtA, tet0PtB, tet0PtC, tet0PtD);

        switch (triI)
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
            case 0:
            {
                n0 = tet0.Sa();
                break;
            }
            case 1:
            {
                n0 = tet0.Sb();
                break;
            }
            case 2:
            {
                n0 = tet0.Sc();
                break;
            }
            case 3:
            {
                n0 = tet0.Sd();
                break;
            }
            default:
            {
                break;
            }
        }
    }

    if (mag(n0) < SMALL)
    {
        // If the old normal is zero (for example in layer addition)
        // then use the current normal;
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
        n0 = n;
    }

    scalar lambdaNumerator = 0;
    scalar lambdaDenominator = 0;

    vector dP = to - from;
    vector dN = n - n0;
    vector dB = b - b0;
    vector dS = from - b0;

    if (mag(dN) > SMALL)
    {
        scalar a = (dP - dB) & dN;
        scalar b = ((dP - dB) & n0) + (dS & dN);
        scalar c = dS & n0;

        if (mag(a) > SMALL)
        {

            // Solve quadratic for lambda
            scalar discriminant = sqr(b) - 4.0*a*c;

            if (discriminant < 0)
            {
                // Imaginary roots only - face not crossed
                return GREAT;
            }
            else
248
            {
249
250
251
252
253
254
255
256
257
258
                // Numerical Recipes in C,
                // Second Edition (1992),
                // Section 5.6.
                // q = -0.5*(b + sgn(b)*sqrt(b^2 - 4ac))
                // x1 = q/a
                // x2 = c/q

                scalar q = -0.5*(b + sign(b)*Foam::sqrt(discriminant));

                if (mag(q) < VSMALL)
259
                {
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
                    // If q is zero, then l1 = q/a is the required
                    // value of lambda, and is zero.

                    return 0.0;
                }

                scalar l1 = q/a;
                scalar l2 = c/q;

                // There will be two roots, a big one and a little
                // one, choose the little one.

                if (mag(l1) < mag(l2))
                {
                    return l1;
                }
                else
                {
                    return l2;
279
                }
280
281
            }
        }
282
283
        {
            // When a is zero, solve the first order polynomial
284

285
286
287
288
289
290
291
292
            lambdaNumerator = -c;
            lambdaDenominator = b;
        }
    }
    else
    {
        // when n = n0 is zero, there is no plane rotation, solve the
        // first order polynomial
293

294
295
        lambdaNumerator = -(dS & n0);
        lambdaDenominator = ((dP - dB) & n0);
296

297
    }
298

299
300
301
    if (mag(lambdaDenominator) < SMALL)
    {
        if (mag(lambdaNumerator) < SMALL)
302
        {
303
304
305
306
307
308
            // Track starts on the face, and is potentially
            // parallel to it.  +-SMALL)/+-SMALL is not a good
            // comparison, return 0.0, in anticipation of tet
            // centre correction.

            return 0.0;
309
310
311
        }
        else
        {
312
313
314
315
316
317
318
319
320
321
322
323
324
            if (mag((to - from)) < SMALL)
            {
                // Zero length track, not along the face, face
                // cannot be crossed.

                return GREAT;
            }
            else
            {
                // Trajectory is non-zero and parallel to face

                lambdaDenominator = sign(lambdaDenominator)*SMALL;
            }
325
        }
326
327
328
329
330
331
332
333
334
335
    }

    return lambdaNumerator/lambdaDenominator;
}


template<class ParticleType>
inline void Foam::Particle<ParticleType>::tetNeighbour(label triI)
{
    const polyMesh& mesh = cloud_.polyMesh_;
336

337
338
    const labelList& pOwner = mesh.faceOwner();
    const faceList& pFaces = mesh.faces();
339

340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
    bool own = (pOwner[tetFaceI_] == cellI_);

    const Foam::face& f = pFaces[tetFaceI_];

    label tetBasePtI = mesh.tetBasePtIs()[tetFaceI_];

    if (tetBasePtI == -1)
    {
        FatalErrorIn
        (
            "inline void Foam::Particle<ParticleType>::tetNeighbour(label triI)"
        )
            << "No base point for face " << tetFaceI_ << ", " << f
            << ", produces a valid tet decomposition."
            << abort(FatalError);
    }

357
358
359
    label facePtI = (tetPtI_ + tetBasePtI) % f.size();
    label otherFacePtI = f.fcIndex(facePtI);

360
361
362
    switch (triI)
    {
        case 0:
363
        {
364
365
366
367
368
            // Crossing this triangle changes tet to that in the
            // neighbour cell over tetFaceI

            // Modification of cellI_ will happen by other indexing,
            // tetFaceI_ and tetPtI don't change.
369

370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
            break;
        }
        case 1:
        {
            crossEdgeConnectedFace
            (
                cellI_,
                tetFaceI_,
                tetPtI_,
                Foam::edge(f[facePtI], f[otherFacePtI])
            );

            break;
        }
        case 2:
        {
            if (own)
            {
                if (tetPtI_ < f.size() - 2)
389
                {
390
                    tetPtI_ = f.fcIndex(tetPtI_);
391
392
393
                }
                else
                {
394
395
396
397
398
399
400
                    crossEdgeConnectedFace
                    (
                        cellI_,
                        tetFaceI_,
                        tetPtI_,
                        Foam::edge(f[tetBasePtI], f[otherFacePtI])
                    );
401
402
403
404
                }
            }
            else
            {
405
406
407
408
409
410
411
412
413
414
415
416
417
418
                if (tetPtI_ > 1)
                {
                    tetPtI_ = f.rcIndex(tetPtI_);
                }
                else
                {
                    crossEdgeConnectedFace
                    (
                        cellI_,
                        tetFaceI_,
                        tetPtI_,
                        Foam::edge(f[tetBasePtI], f[facePtI])
                    );
                }
419
            }
420
421

            break;
422
        }
423
        case 3:
424
        {
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
            if (own)
            {
                if (tetPtI_ > 1)
                {
                    tetPtI_ = f.rcIndex(tetPtI_);
                }
                else
                {
                    crossEdgeConnectedFace
                    (
                        cellI_,
                        tetFaceI_,
                        tetPtI_,
                        Foam::edge(f[tetBasePtI], f[facePtI])
                    );
                }
            }
            else
443
            {
444
                if (tetPtI_ < f.size() - 2)
445
                {
446
                    tetPtI_ = f.fcIndex(tetPtI_);
447
448
449
                }
                else
                {
450
451
452
453
454
455
456
                    crossEdgeConnectedFace
                    (
                        cellI_,
                        tetFaceI_,
                        tetPtI_,
                        Foam::edge(f[tetBasePtI], f[otherFacePtI])
                    );
457
458
459
                }
            }

460
461
462
463
464
465
466
467
468
469
470
471
472
473
            break;
        }
        default:
        {
            FatalErrorIn
            (
                "inline void "
                "Foam::Particle<ParticleType>::tetNeighbour(label triI)"
            )
                << "Tet tri face index error, can only be 0..3, supplied "
                << triI << nl
                << abort(FatalError);

            break;
474
475
476
477
478
        }
    }
}

template<class ParticleType>
479
inline void Foam::Particle<ParticleType>::crossEdgeConnectedFace
480
(
481
482
483
484
485
    const label& cellI,
    label& tetFaceI,
    label& tetPtI,
    const edge& e
)
486
487
488
{
    const polyMesh& mesh = cloud_.polyMesh_;

489
490
491
492
493
494
    const faceList& pFaces = mesh.faces();
    const cellList& pCells = mesh.cells();

    const Foam::face& f = pFaces[tetFaceI];

    const Foam::cell& thisCell = pCells[cellI];
495

496
    forAll(thisCell, cFI)
497
    {
498
499
        // Loop over all other faces of this cell and
        // find the one that shares this edge
500

501
502
503
        label fI = thisCell[cFI];

        if (tetFaceI == fI)
504
        {
505
            continue;
506
507
        }

508
        const Foam::face& otherFace = pFaces[fI];
509

510
511
512
        label edDir = otherFace.edgeDirection(e);

        if (edDir == 0)
513
        {
514
            continue;
515
516
517
        }
        else
        {
518
519
            //Found edge on other face
            tetFaceI = fI;
520

521
            label eIndex = -1;
522

523
524
525
526
            if (edDir == 1)
            {
                // Edge is in the forward circulation of this face, so
                // work with the start point of the edge
527

528
529
530
531
532
533
534
                eIndex = findIndex(otherFace, e.start());
            }
            else
            {
                // edDir == -1, so the edge is in the reverse
                // circulation of this face, so work with the end
                // point of the edge
535

536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
                eIndex = findIndex(otherFace, e.end());
            }

            label tetBasePtI = mesh.tetBasePtIs()[fI];

            if (tetBasePtI == -1)
            {
                FatalErrorIn
                (
                    "inline void "
                    "Foam::Particle<ParticleType>::crossEdgeConnectedFace"
                    "("
                        "const label& cellI,"
                        "label& tetFaceI,"
                        "label& tetPtI,"
                        "const edge& e"
                    ")"
                )
                    << "No base point for face " << fI << ", " << f
                    << ", produces a decomposition that has a minimum "
                    << "volume greater than tolerance."
                    << abort(FatalError);
            }

            // Find eIndex relative to the base point on new face
            eIndex -= tetBasePtI;

            if (neg(eIndex))
            {
                eIndex = (eIndex + otherFace.size()) % otherFace.size();
            }

            if (eIndex == 0)
            {
                // The point is the base point, so this is first tet
                // in the face circulation

                tetPtI = 1;
            }
            else if (eIndex == otherFace.size() - 1)
            {
                // The point is the last before the base point, so
                // this is the last tet in the face circulation

                tetPtI = otherFace.size() - 2;
            }
            else
            {
                tetPtI = eIndex;
            }

            break;
        }
    }
590
591
592
593
}


template<class ParticleType>
594
inline void Foam::Particle<ParticleType>::hitWallFaces
595
(
596
597
598
599
600
    const vector& from,
    const vector& to,
    scalar& lambdaMin,
    tetIndices& closestTetIs
)
601
{
602
603
604
605
606
607
608
609
610
    if (!(cloud_.hasWallImpactDistance() && cloud_.cellHasWallFaces()[cellI_]))
    {
        return;
    }

    const polyMesh& mesh = cloud_.polyMesh_;
    const faceList& pFaces = mesh.faces();

    const ParticleType& p = static_cast<const ParticleType&>(*this);
611

612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
    const Foam::cell& thisCell = mesh.cells()[cellI_];

    const polyBoundaryMesh& patches = mesh.boundaryMesh();

    forAll(thisCell, cFI)
    {
        label fI = thisCell[cFI];

        if (cloud_.internalFace(fI))
        {
            continue;
        }

        label patchI = patches.patchID()[fI - mesh.nInternalFaces()];

        if (isA<wallPolyPatch>(patches[patchI]))
        {
            // Get the decomposition of this wall face

            const List<tetIndices>& faceTetIs =
                polyMeshTetDecomposition::faceTetIndices(mesh, fI, cellI_);

            const Foam::face& f = pFaces[fI];

            forAll(faceTetIs, tI)
            {
                const tetIndices& tetIs = faceTetIs[tI];

                triPointRef tri = tetIs.faceTri(mesh);

                vector n = tri.normal();

                vector nHat = n/mag(n);

                // Radius of particle with respect to this wall face
                // triangle.  Assuming that the wallImpactDistance
                // does not change as the particle or the mesh moves
                // forward in time.
                scalar r = p.wallImpactDistance(nHat);

                vector toPlusRNHat = to + r*nHat;

                // triI = 0 because it is the cell face tri of the tet
                // we are concerned with.
                scalar tetClambda = tetLambda
                (
                    tetIs.tet(mesh).centre(),
                    toPlusRNHat,
                    0,
                    n,
                    f[tetIs.faceBasePt()],
                    cellI_,
                    fI,
                    tetIs.tetPt()
                );

                if ((tetClambda <= 0.0) || (tetClambda >= 1.0))
                {
                    // toPlusRNHat is not on the outside of the plane of
                    // the wall face tri, the tri cannot be hit.
                    continue;
                }

                // Check if the actual trajectory of the near-tri
                // points intersects the triangle.

                vector fromPlusRNHat = from + r*nHat;

                // triI = 0 because it is the cell face tri of the tet
                // we are concerned with.
                scalar lambda = tetLambda
                (
                    fromPlusRNHat,
                    toPlusRNHat,
                    0,
                    n,
                    f[tetIs.faceBasePt()],
                    cellI_,
                    fI,
                    tetIs.tetPt()
                );

                pointHit hitInfo(vector::zero);

                if (mesh.moving())
                {
                    // For a moving mesh, the position of wall
                    // triangle needs to be moved in time to be
                    // consistent with the moment defined by the
                    // current value of stepFraction and the value of
                    // lambda just calculated.

                    // Total fraction thought the timestep of the
                    // motion, including stepFraction before the
                    // current tracking step and the current
                    // lambda
                    // i.e.
                    // let s = stepFraction, l = lambda
                    // Motion of x in time:
                    // |-----------------|---------|---------|
                    // x00               x0        xi        x
                    //
                    // where xi is the correct value of x at the required
                    // tracking instant.
                    //
                    // x0 = x00 + s*(x - x00) = s*x + (1 - s)*x00
                    //
                    // i.e. the motion covered by previous tracking portions
                    // within this timestep, and
                    //
                    // xi = x0 + l*(x - x0)
                    //    = l*x + (1 - l)*x0
                    //    = l*x + (1 - l)*(s*x + (1 - s)*x00)
                    //    = (s + l - s*l)*x + (1 - (s + l - s*l))*x00
                    //
                    // let m = (s + l - s*l)
                    //
                    // xi = m*x + (1 - m)*x00 = x00 + m*(x - x00);
                    //
                    // In the same form as before.

                    // Clip lambda to 0.0-1.0 to ensure that sensible
                    // positions are used for triangle intersections.
                    scalar lam = max(0.0, min(1.0, lambda));

                    scalar m = stepFraction_ + lam - (stepFraction_*lam);

                    triPointRef tri00 = tetIs.oldFaceTri(mesh);

                    // Use SMALL positive tolerance to make the triangle
                    // slightly "fat" to improve robustness.  Intersection
                    // is calculated as the ray (from + r*nHat) -> (to +
                    // r*nHat).

                    point tPtA = tri00.a() + m*(tri.a() - tri00.a());
                    point tPtB = tri00.b() + m*(tri.b() - tri00.b());
                    point tPtC = tri00.c() + m*(tri.c() - tri00.c());

                    triPointRef t(tPtA, tPtB, tPtC);

                    // The point fromPlusRNHat + m*(to - from) is on the
                    // plane of the triangle.  Determine the
                    // intersection with this triangle by testing if
                    // this point is inside or outside of the triangle.
                    hitInfo = t.intersection
                    (
                        fromPlusRNHat + m*(to - from),
                        t.normal(),
                        intersection::FULL_RAY,
                        SMALL
                    );
                }
                else
                {
                    // Use SMALL positive tolerance to make the triangle
                    // slightly "fat" to improve robustness.  Intersection
                    // is calculated as the ray (from + r*nHat) -> (to +
                    // r*nHat).
                    hitInfo = tri.intersection
                    (
                        fromPlusRNHat,
                        (to - from),
                        intersection::FULL_RAY,
                        SMALL
                    );
                }

                if (hitInfo.hit())
                {
                    if (lambda < lambdaMin)
                    {
                        lambdaMin = lambda;

                        faceI_ = fI;

                        closestTetIs = tetIs;
                    }
                }
            }
        }
    }
793
794
795
}


796
// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
797
798

template<class ParticleType>
799
inline Foam::Particle<ParticleType>::trackData::trackData
800
801
802
803
804
805
806
(
    Cloud<ParticleType>& cloud
)
:
    cloud_(cloud)
{}

807

808
809
// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

810
template<class ParticleType>
811
812
inline Foam::Cloud<ParticleType>&
Foam::Particle<ParticleType>::trackData::cloud()
813
814
815
816
817
818
{
    return cloud_;
}


template<class ParticleType>
819
820
inline const Foam::Cloud<ParticleType>&
Foam::Particle<ParticleType>::cloud() const
821
822
823
824
825
826
{
    return cloud_;
}


template<class ParticleType>
827
inline const Foam::vector& Foam::Particle<ParticleType>::position() const
828
829
830
831
832
833
{
    return position_;
}


template<class ParticleType>
834
inline Foam::vector& Foam::Particle<ParticleType>::position()
835
836
837
838
839
840
{
    return position_;
}


template<class ParticleType>
841
inline Foam::label Foam::Particle<ParticleType>::cell() const
842
{
843
    return cellI_;
844
845
}

846

847
template<class ParticleType>
848
inline Foam::label& Foam::Particle<ParticleType>::cell()
849
{
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
    return cellI_;
}


template<class ParticleType>
inline Foam::label Foam::Particle<ParticleType>::tetFace() const
{
    return tetFaceI_;
}


template<class ParticleType>
inline Foam::label& Foam::Particle<ParticleType>::tetFace()
{
    return tetFaceI_;
}


template<class ParticleType>
inline Foam::label Foam::Particle<ParticleType>::tetPt() const
{
    return tetPtI_;
}


template<class ParticleType>
inline Foam::label& Foam::Particle<ParticleType>::tetPt()
{
    return tetPtI_;
}


template<class ParticleType>
inline Foam::tetIndices Foam::Particle<ParticleType>::currentTetIndices() const
{
    return tetIndices(cellI_, tetFaceI_, tetPtI_, cloud_.polyMesh_);
}


template<class ParticleType>
inline Foam::tetPointRef Foam::Particle<ParticleType>::currentTet() const
{
    return currentTetIndices().tet(cloud_.polyMesh_);
}


template<class ParticleType>
inline Foam::vector Foam::Particle<ParticleType>::normal() const
{
    return currentTetIndices().faceTri(cloud_.polyMesh_).normal();
}


template<class ParticleType>
inline Foam::vector
Foam::Particle<ParticleType>::oldNormal() const
{
    return currentTetIndices().oldFaceTri(cloud_.polyMesh_).normal();
908
909
}

910
911

template<class ParticleType>
912
inline Foam::label Foam::Particle<ParticleType>::face() const
913
{
914
    return faceI_;
915
916
917
}


918
919
920
template<class ParticleType>
inline Foam::label& Foam::Particle<ParticleType>::face()
{
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
    return faceI_;
}


template<class ParticleType>
inline void Foam::Particle<ParticleType>::initCellFacePt()
{
    if (cellI_ == -1)
    {
        cloud_.findCellFacePt
        (
            position_,
            cellI_,
            tetFaceI_,
            tetPtI_
        );

        if (cellI_ == -1)
        {
            FatalErrorIn
            (
                "Foam::Particle<ParticleType>::Particle"
                "("
                    "void Foam::Particle<ParticleType>::initCellFacePt()"
                ")"
            )   << "cell, tetFace and tetPt search failure at position "
                << position_ << nl
                << abort(FatalError);
        }
    }
    else
    {
        cloud_.findFacePt(cellI_, position_, tetFaceI_, tetPtI_);

        if (tetFaceI_ == -1 || tetPtI_ == -1)
        {
            label oldCellI = cellI_;

            cloud_.findCellFacePt
            (
                position_,
                cellI_,
                tetFaceI_,
                tetPtI_
            );

            if (cellI_ == -1 || tetFaceI_ == -1 || tetPtI_ == -1)
            {
                // The particle has entered this function with a cell
                // number, but hasn't been able to find a cell to
                // occupy.

973
                if(!cloud_.polyMesh_.pointInCellBB(position_, oldCellI, 0.1))
974
                {
975
976
977
978
                    // If the position is not inside the (slightly
                    // extended) bound-box of the cell that it thought
                    // it should be in, then this is considered an
                    // error.
979
980
981
982
983
984
985
986
987
988

                    FatalErrorIn
                    (
                        "Foam::Particle<ParticleType>::Particle"
                        "("
                            "void "
                            "Foam::Particle<ParticleType>::initCellFacePt()"
                        ")"
                    )   << "cell, tetFace and tetPt search failure at position "
                        << position_ << nl
989
                        << "for requested cell " << oldCellI << nl
990
991
992
                        << abort(FatalError);
                }

993
994
995
996
997
998
999
1000
1001
1002
                // The position is in the (slightly extended)
                // bound-box of the cell.  This situation may arise
                // because the face decomposition of the cell is not
                // the same as when the particle acquired the cell
                // index.  For example, it has been read into a mesh
                // that has made a different face base-point decision
                // for a boundary face and now this particle is in a
                // position that is not in the mesh.  Gradually move
                // the particle towards the centre of the cell that it
                // thought that it was in.
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088

                cellI_ = oldCellI;

                point newPosition = position_;

                const point& cC = cloud_.polyMesh_.cellCentres()[cellI_];

                label trap(1.0/Cloud<ParticleType>::trackingCorrectionTol + 1);

                label iterNo = 0;

                do
                {
                    newPosition +=
                        Cloud<ParticleType>::trackingCorrectionTol
                       *(cC - position_);

                    cloud_.findFacePt
                    (
                        cellI_,
                        newPosition,
                        tetFaceI_,
                        tetPtI_
                    );

                    iterNo++;

                } while (tetFaceI_ < 0  && iterNo <= trap);

                if(tetFaceI_ == -1)
                {
                    FatalErrorIn
                    (
                        "Foam::Particle<ParticleType>::Particle"
                        "("
                        "void "
                        "Foam::Particle<ParticleType>::initCellFacePt()"
                        ")"
                    )   << "cell, tetFace and tetPt search failure at position "
                        << position_ << nl
                        << abort(FatalError);
                }

                WarningIn
                (
                    "Foam::Particle<ParticleType>::Particle"
                    "("
                        "void Foam::Particle<ParticleType>::initCellFacePt()"
                    ")"
                )
                    << "Particle moved from " << position_
                    << " to " << newPosition
                    << " in cell " << cellI_
                    << " tetFace " << tetFaceI_
                    << " tetPt " << tetPtI_ << nl
                    << "    (A fraction of "
                    << 1.0 - mag(cC - newPosition)/mag(cC - position_)
                    << " of the distance to the cell centre)"
                    << " because a decomposition tetFace and tetPt "
                    << "could not be found."
                    << endl;

                position_ = newPosition;
            }

            if (cellI_ != oldCellI)
            {
                WarningIn
                (
                    "Foam::Particle<ParticleType>::Particle"
                    "("
                        "void Foam::Particle<ParticleType>::initCellFacePt()"
                    ")"
                )
                    << "Particle at position " << position_
                    << " searched for a cell, tetFace and tetPt." << nl
                    << "    Found"
                    << " cell " << cellI_
                    << " tetFace " << tetFaceI_
                    << " tetPt " << tetPtI_ << nl
                    << "    This is a different cell to that which was supplied"
                    << " (" << oldCellI << ")." << nl
                    << endl;
            }
        }
    }
1089
1090
1091
}


1092
template<class ParticleType>
1093
inline bool Foam::Particle<ParticleType>::onBoundary() const
1094
{
1095
    return faceI_ != -1 && faceI_ >= cloud_.pMesh().nInternalFaces();
1096
1097
1098
1099
}


template<class ParticleType>
1100
inline Foam::scalar& Foam::Particle<ParticleType>::stepFraction()
1101
1102
1103
1104
1105
1106
{
    return stepFraction_;
}


template<class ParticleType>
1107
inline Foam::scalar Foam::Particle<ParticleType>::stepFraction() const
1108
1109
1110
1111
1112
1113
{
    return stepFraction_;
}


template<class ParticleType>
1114
1115
1116
1117
1118
1119
inline Foam::label Foam::Particle<ParticleType>::origProc() const
{
    return origProc_;
}


1120
1121
1122
1123
1124
1125
1126
template<class ParticleType>
inline Foam::label& Foam::Particle<ParticleType>::origProc()
{
    return origProc_;
}


1127
1128
1129
1130
1131
1132
template<class ParticleType>
inline Foam::label Foam::Particle<ParticleType>::origId() const
{
    return origId_;
}

1133
1134
1135
1136
1137
1138
1139

template<class ParticleType>
inline Foam::label& Foam::Particle<ParticleType>::origId()
{
    return origId_;
}

1140
1141
1142

template<class ParticleType>
inline bool Foam::Particle<ParticleType>::softImpact() const
1143
1144
1145
1146
1147
{
    return false;
}


1148
template<class ParticleType>
1149
inline Foam::scalar Foam::Particle<ParticleType>::currentTime() const
1150
1151
1152
{
    return
        cloud_.pMesh().time().value()
graham's avatar
graham committed
1153
      + stepFraction_*cloud_.pMesh().time().deltaTValue();
1154
1155
1156
}


1157
template<class ParticleType>
1158
inline Foam::label Foam::Particle<ParticleType>::patch(const label faceI) const
1159
{
1160
    return cloud_.facePatch(faceI);
1161
1162
1163
1164
}


template<class ParticleType>
1165
inline Foam::label Foam::Particle<ParticleType>::patchFace
1166
(
1167
1168
    const label patchI,
    const label faceI
1169
1170
) const
{
1171
    return cloud_.patchFace(patchI, faceI);
1172
1173
1174
1175
}


template<class ParticleType>
1176
1177
inline Foam::scalar
Foam::Particle<ParticleType>::wallImpactDistance(const vector&) const
1178
1179
1180
1181
1182
1183
{
    return 0.0;
}


template<class ParticleType>
1184
inline Foam::label Foam::Particle<ParticleType>::faceInterpolation() const
1185
{
1186
    return faceI_;
1187
1188
1189
1190
}


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